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.
45 #include "gromacs/fileio/confio.h"
46 #include "gromacs/fileio/mtxio.h"
47 #include "gromacs/fileio/trajectory_writing.h"
48 #include "gromacs/imd/imd.h"
49 #include "gromacs/legacyheaders/bonded-threading.h"
50 #include "gromacs/legacyheaders/constr.h"
51 #include "gromacs/legacyheaders/domdec.h"
52 #include "gromacs/legacyheaders/force.h"
53 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
54 #include "gromacs/legacyheaders/macros.h"
55 #include "gromacs/legacyheaders/md_logging.h"
56 #include "gromacs/legacyheaders/md_support.h"
57 #include "gromacs/legacyheaders/mdatoms.h"
58 #include "gromacs/legacyheaders/mdebin.h"
59 #include "gromacs/legacyheaders/mdrun.h"
60 #include "gromacs/legacyheaders/names.h"
61 #include "gromacs/legacyheaders/network.h"
62 #include "gromacs/legacyheaders/nrnb.h"
63 #include "gromacs/legacyheaders/ns.h"
64 #include "gromacs/legacyheaders/pme.h"
65 #include "gromacs/legacyheaders/sim_util.h"
66 #include "gromacs/legacyheaders/tgroup.h"
67 #include "gromacs/legacyheaders/txtdump.h"
68 #include "gromacs/legacyheaders/typedefs.h"
69 #include "gromacs/legacyheaders/update.h"
70 #include "gromacs/legacyheaders/vsite.h"
71 #include "gromacs/legacyheaders/types/commrec.h"
72 #include "gromacs/linearalgebra/sparsematrix.h"
73 #include "gromacs/math/vec.h"
74 #include "gromacs/pbcutil/mshift.h"
75 #include "gromacs/pbcutil/pbc.h"
76 #include "gromacs/timing/wallcycle.h"
77 #include "gromacs/timing/walltime_accounting.h"
78 #include "gromacs/topology/mtop_util.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,
314 gmx_wallcycle_t wcycle)
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, 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, wcycle);
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, 1.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(inputrec, fr, 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, 1.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 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
794 m_add(force_vir, shake_vir, vir);
795 wallcycle_stop(wcycle, ewcCONSTR);
799 copy_mat(force_vir, vir);
803 enerd->term[F_PRES] =
804 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
806 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
808 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
810 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
814 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
816 em_state_t *s_min, em_state_t *s_b)
820 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
822 unsigned char *grpnrFREEZE;
826 fprintf(debug, "Doing reorder_partsum\n");
832 cgs_gl = dd_charge_groups_global(cr->dd);
833 index = cgs_gl->index;
835 /* Collect fm in a global vector fmg.
836 * This conflicts with the spirit of domain decomposition,
837 * but to fully optimize this a much more complicated algorithm is required.
839 snew(fmg, mtop->natoms);
841 ncg = s_min->s.ncg_gl;
842 cg_gl = s_min->s.cg_gl;
844 for (c = 0; c < ncg; c++)
849 for (a = a0; a < a1; a++)
851 copy_rvec(fm[i], fmg[a]);
855 gmx_sum(mtop->natoms*3, fmg[0], cr);
857 /* Now we will determine the part of the sum for the cgs in state s_b */
859 cg_gl = s_b->s.cg_gl;
863 grpnrFREEZE = mtop->groups.grpnr[egcFREEZE];
864 for (c = 0; c < ncg; c++)
869 for (a = a0; a < a1; a++)
871 if (mdatoms->cFREEZE && grpnrFREEZE)
875 for (m = 0; m < DIM; m++)
877 if (!opts->nFreeze[gf][m])
879 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
891 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
893 em_state_t *s_min, em_state_t *s_b)
899 /* This is just the classical Polak-Ribiere calculation of beta;
900 * it looks a bit complicated since we take freeze groups into account,
901 * and might have to sum it in parallel runs.
904 if (!DOMAINDECOMP(cr) ||
905 (s_min->s.ddp_count == cr->dd->ddp_count &&
906 s_b->s.ddp_count == cr->dd->ddp_count))
912 /* This part of code can be incorrect with DD,
913 * since the atom ordering in s_b and s_min might differ.
915 for (i = 0; i < mdatoms->homenr; i++)
917 if (mdatoms->cFREEZE)
919 gf = mdatoms->cFREEZE[i];
921 for (m = 0; m < DIM; m++)
923 if (!opts->nFreeze[gf][m])
925 sum += (fb[i][m] - fm[i][m])*fb[i][m];
932 /* We need to reorder cgs while summing */
933 sum = reorder_partsum(cr, opts, mdatoms, mtop, s_min, s_b);
937 gmx_sumd(1, &sum, cr);
940 return sum/sqr(s_min->fnorm);
943 double do_cg(FILE *fplog, t_commrec *cr,
944 int nfile, const t_filenm fnm[],
945 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
946 int gmx_unused nstglobalcomm,
947 gmx_vsite_t *vsite, gmx_constr_t constr,
948 int gmx_unused stepout,
949 t_inputrec *inputrec,
950 gmx_mtop_t *top_global, t_fcdata *fcd,
951 t_state *state_global,
953 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
954 gmx_edsam_t gmx_unused ed,
956 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
957 gmx_membed_t gmx_unused membed,
958 real gmx_unused cpt_period, real gmx_unused max_hours,
959 const char gmx_unused *deviceOptions,
961 unsigned long gmx_unused Flags,
962 gmx_walltime_accounting_t walltime_accounting)
964 const char *CG = "Polak-Ribiere Conjugate Gradients";
966 em_state_t *s_min, *s_a, *s_b, *s_c;
968 gmx_enerdata_t *enerd;
970 gmx_global_stat_t gstat;
972 rvec *f_global, *p, *sf, *sfm;
973 double gpa, gpb, gpc, tmp, sum[2], minstep;
976 real a, b, c, beta = 0.0;
980 gmx_bool converged, foundlower;
982 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
984 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
986 int i, m, gf, step, nminstep;
991 s_min = init_em_state();
992 s_a = init_em_state();
993 s_b = init_em_state();
994 s_c = init_em_state();
996 /* Init em and store the local state in s_min */
997 init_em(fplog, CG, cr, inputrec,
998 state_global, top_global, s_min, &top, &f, &f_global,
999 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1000 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
1002 /* Print to log file */
1003 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1005 /* Max number of steps */
1006 number_steps = inputrec->nsteps;
1010 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1014 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1017 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1018 /* do_force always puts the charge groups in the box and shifts again
1019 * We do not unshift, so molecules are always whole in congrad.c
1021 evaluate_energy(fplog, cr,
1022 top_global, s_min, top,
1023 inputrec, nrnb, wcycle, gstat,
1024 vsite, constr, fcd, graph, mdatoms, fr,
1025 mu_tot, enerd, vir, pres, -1, TRUE);
1030 /* Copy stuff to the energy bin for easy printing etc. */
1031 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1032 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1033 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1035 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1036 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1037 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1041 /* Estimate/guess the initial stepsize */
1042 stepsize = inputrec->em_stepsize/s_min->fnorm;
1046 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1047 s_min->fmax, s_min->a_fmax+1);
1048 fprintf(stderr, " F-Norm = %12.5e\n",
1049 s_min->fnorm/sqrt(state_global->natoms));
1050 fprintf(stderr, "\n");
1051 /* and copy to the log file too... */
1052 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1053 s_min->fmax, s_min->a_fmax+1);
1054 fprintf(fplog, " F-Norm = %12.5e\n",
1055 s_min->fnorm/sqrt(state_global->natoms));
1056 fprintf(fplog, "\n");
1058 /* Start the loop over CG steps.
1059 * Each successful step is counted, and we continue until
1060 * we either converge or reach the max number of steps.
1063 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1066 /* start taking steps in a new direction
1067 * First time we enter the routine, beta=0, and the direction is
1068 * simply the negative gradient.
1071 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1076 for (i = 0; i < mdatoms->homenr; i++)
1078 if (mdatoms->cFREEZE)
1080 gf = mdatoms->cFREEZE[i];
1082 for (m = 0; m < DIM; m++)
1084 if (!inputrec->opts.nFreeze[gf][m])
1086 p[i][m] = sf[i][m] + beta*p[i][m];
1087 gpa -= p[i][m]*sf[i][m];
1088 /* f is negative gradient, thus the sign */
1097 /* Sum the gradient along the line across CPUs */
1100 gmx_sumd(1, &gpa, cr);
1103 /* Calculate the norm of the search vector */
1104 get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1106 /* Just in case stepsize reaches zero due to numerical precision... */
1109 stepsize = inputrec->em_stepsize/pnorm;
1113 * Double check the value of the derivative in the search direction.
1114 * If it is positive it must be due to the old information in the
1115 * CG formula, so just remove that and start over with beta=0.
1116 * This corresponds to a steepest descent step.
1121 step--; /* Don't count this step since we are restarting */
1122 continue; /* Go back to the beginning of the big for-loop */
1125 /* Calculate minimum allowed stepsize, before the average (norm)
1126 * relative change in coordinate is smaller than precision
1129 for (i = 0; i < mdatoms->homenr; i++)
1131 for (m = 0; m < DIM; m++)
1133 tmp = fabs(s_min->s.x[i][m]);
1142 /* Add up from all CPUs */
1145 gmx_sumd(1, &minstep, cr);
1148 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1150 if (stepsize < minstep)
1156 /* Write coordinates if necessary */
1157 do_x = do_per_step(step, inputrec->nstxout);
1158 do_f = do_per_step(step, inputrec->nstfout);
1160 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1161 top_global, inputrec, step,
1162 s_min, state_global, f_global);
1164 /* Take a step downhill.
1165 * In theory, we should minimize the function along this direction.
1166 * That is quite possible, but it turns out to take 5-10 function evaluations
1167 * for each line. However, we dont really need to find the exact minimum -
1168 * it is much better to start a new CG step in a modified direction as soon
1169 * as we are close to it. This will save a lot of energy evaluations.
1171 * In practice, we just try to take a single step.
1172 * If it worked (i.e. lowered the energy), we increase the stepsize but
1173 * the continue straight to the next CG step without trying to find any minimum.
1174 * If it didn't work (higher energy), there must be a minimum somewhere between
1175 * the old position and the new one.
1177 * Due to the finite numerical accuracy, it turns out that it is a good idea
1178 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1179 * This leads to lower final energies in the tests I've done. / Erik
1181 s_a->epot = s_min->epot;
1183 c = a + stepsize; /* reference position along line is zero */
1185 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1187 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1188 s_min, top, mdatoms, fr, vsite, constr,
1192 /* Take a trial step (new coords in s_c) */
1193 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1194 constr, top, nrnb, wcycle, -1);
1197 /* Calculate energy for the trial step */
1198 evaluate_energy(fplog, cr,
1199 top_global, s_c, top,
1200 inputrec, nrnb, wcycle, gstat,
1201 vsite, constr, fcd, graph, mdatoms, fr,
1202 mu_tot, enerd, vir, pres, -1, FALSE);
1204 /* Calc derivative along line */
1208 for (i = 0; i < mdatoms->homenr; i++)
1210 for (m = 0; m < DIM; m++)
1212 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1215 /* Sum the gradient along the line across CPUs */
1218 gmx_sumd(1, &gpc, cr);
1221 /* This is the max amount of increase in energy we tolerate */
1222 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1224 /* Accept the step if the energy is lower, or if it is not significantly higher
1225 * and the line derivative is still negative.
1227 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1230 /* Great, we found a better energy. Increase step for next iteration
1231 * if we are still going down, decrease it otherwise
1235 stepsize *= 1.618034; /* The golden section */
1239 stepsize *= 0.618034; /* 1/golden section */
1244 /* New energy is the same or higher. We will have to do some work
1245 * to find a smaller value in the interval. Take smaller step next time!
1248 stepsize *= 0.618034;
1254 /* OK, if we didn't find a lower value we will have to locate one now - there must
1255 * be one in the interval [a=0,c].
1256 * The same thing is valid here, though: Don't spend dozens of iterations to find
1257 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1258 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1260 * I also have a safeguard for potentially really patological functions so we never
1261 * take more than 20 steps before we give up ...
1263 * If we already found a lower value we just skip this step and continue to the update.
1271 /* Select a new trial point.
1272 * If the derivatives at points a & c have different sign we interpolate to zero,
1273 * otherwise just do a bisection.
1275 if (gpa < 0 && gpc > 0)
1277 b = a + gpa*(a-c)/(gpc-gpa);
1284 /* safeguard if interpolation close to machine accuracy causes errors:
1285 * never go outside the interval
1287 if (b <= a || b >= c)
1292 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1294 /* Reload the old state */
1295 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1296 s_min, top, mdatoms, fr, vsite, constr,
1300 /* Take a trial step to this new point - new coords in s_b */
1301 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1302 constr, top, nrnb, wcycle, -1);
1305 /* Calculate energy for the trial step */
1306 evaluate_energy(fplog, cr,
1307 top_global, s_b, top,
1308 inputrec, nrnb, wcycle, gstat,
1309 vsite, constr, fcd, graph, mdatoms, fr,
1310 mu_tot, enerd, vir, pres, -1, FALSE);
1312 /* p does not change within a step, but since the domain decomposition
1313 * might change, we have to use cg_p of s_b here.
1318 for (i = 0; i < mdatoms->homenr; i++)
1320 for (m = 0; m < DIM; m++)
1322 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1325 /* Sum the gradient along the line across CPUs */
1328 gmx_sumd(1, &gpb, cr);
1333 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1334 s_a->epot, s_b->epot, s_c->epot, gpb);
1337 epot_repl = s_b->epot;
1339 /* Keep one of the intervals based on the value of the derivative at the new point */
1342 /* Replace c endpoint with b */
1343 swap_em_state(s_b, s_c);
1349 /* Replace a endpoint with b */
1350 swap_em_state(s_b, s_a);
1356 * Stop search as soon as we find a value smaller than the endpoints.
1357 * Never run more than 20 steps, no matter what.
1361 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1364 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1367 /* OK. We couldn't find a significantly lower energy.
1368 * If beta==0 this was steepest descent, and then we give up.
1369 * If not, set beta=0 and restart with steepest descent before quitting.
1379 /* Reset memory before giving up */
1385 /* Select min energy state of A & C, put the best in B.
1387 if (s_c->epot < s_a->epot)
1391 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1392 s_c->epot, s_a->epot);
1394 swap_em_state(s_b, s_c);
1402 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1403 s_a->epot, s_c->epot);
1405 swap_em_state(s_b, s_a);
1415 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1418 swap_em_state(s_b, s_c);
1423 /* new search direction */
1424 /* beta = 0 means forget all memory and restart with steepest descents. */
1425 if (nstcg && ((step % nstcg) == 0))
1431 /* s_min->fnorm cannot be zero, because then we would have converged
1435 /* Polak-Ribiere update.
1436 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1438 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1440 /* Limit beta to prevent oscillations */
1441 if (fabs(beta) > 5.0)
1447 /* update positions */
1448 swap_em_state(s_min, s_b);
1451 /* Print it if necessary */
1456 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1457 step, s_min->epot, s_min->fnorm/sqrt(state_global->natoms),
1458 s_min->fmax, s_min->a_fmax+1);
1460 /* Store the new (lower) energies */
1461 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1462 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1463 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1465 do_log = do_per_step(step, inputrec->nstlog);
1466 do_ene = do_per_step(step, inputrec->nstenergy);
1468 /* Prepare IMD energy record, if bIMD is TRUE. */
1469 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1473 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1475 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1476 do_log ? fplog : NULL, step, step, eprNORMAL,
1477 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1480 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1481 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
1483 IMD_send_positions(inputrec->imd);
1486 /* Stop when the maximum force lies below tolerance.
1487 * If we have reached machine precision, converged is already set to true.
1489 converged = converged || (s_min->fmax < inputrec->em_tol);
1491 } /* End of the loop */
1493 /* IMD cleanup, if bIMD is TRUE. */
1494 IMD_finalize(inputrec->bIMD, inputrec->imd);
1498 step--; /* we never took that last step in this case */
1501 if (s_min->fmax > inputrec->em_tol)
1505 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1506 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1513 /* If we printed energy and/or logfile last step (which was the last step)
1514 * we don't have to do it again, but otherwise print the final values.
1518 /* Write final value to log since we didn't do anything the last step */
1519 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1521 if (!do_ene || !do_log)
1523 /* Write final energy file entries */
1524 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1525 !do_log ? fplog : NULL, step, step, eprNORMAL,
1526 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1530 /* Print some stuff... */
1533 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1537 * For accurate normal mode calculation it is imperative that we
1538 * store the last conformation into the full precision binary trajectory.
1540 * However, we should only do it if we did NOT already write this step
1541 * above (which we did if do_x or do_f was true).
1543 do_x = !do_per_step(step, inputrec->nstxout);
1544 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1546 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1547 top_global, inputrec, step,
1548 s_min, state_global, f_global);
1550 fnormn = s_min->fnorm/sqrt(state_global->natoms);
1554 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1555 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1556 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1557 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1559 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1562 finish_em(cr, outf, walltime_accounting, wcycle);
1564 /* To print the actual number of steps we needed somewhere */
1565 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1568 } /* That's all folks */
1571 double do_lbfgs(FILE *fplog, t_commrec *cr,
1572 int nfile, const t_filenm fnm[],
1573 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
1574 int gmx_unused nstglobalcomm,
1575 gmx_vsite_t *vsite, gmx_constr_t constr,
1576 int gmx_unused stepout,
1577 t_inputrec *inputrec,
1578 gmx_mtop_t *top_global, t_fcdata *fcd,
1581 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1582 gmx_edsam_t gmx_unused ed,
1584 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1585 gmx_membed_t gmx_unused membed,
1586 real gmx_unused cpt_period, real gmx_unused max_hours,
1587 const char gmx_unused *deviceOptions,
1589 unsigned long gmx_unused Flags,
1590 gmx_walltime_accounting_t walltime_accounting)
1592 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1594 gmx_localtop_t *top;
1595 gmx_enerdata_t *enerd;
1597 gmx_global_stat_t gstat;
1600 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1601 double stepsize, gpa, gpb, gpc, tmp, minstep;
1602 real *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1603 real *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1604 real a, b, c, maxdelta, delta;
1605 real diag, Epot0, Epot, EpotA, EpotB, EpotC;
1606 real dgdx, dgdg, sq, yr, beta;
1608 gmx_bool converged, first;
1611 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1613 int start, end, number_steps;
1615 int i, k, m, n, nfmax, gf, step;
1622 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1627 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).");
1630 n = 3*state->natoms;
1631 nmaxcorr = inputrec->nbfgscorr;
1633 /* Allocate memory */
1634 /* Use pointers to real so we dont have to loop over both atoms and
1635 * dimensions all the time...
1636 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1637 * that point to the same memory.
1650 snew(rho, nmaxcorr);
1651 snew(alpha, nmaxcorr);
1654 for (i = 0; i < nmaxcorr; i++)
1660 for (i = 0; i < nmaxcorr; i++)
1669 init_em(fplog, LBFGS, cr, inputrec,
1670 state, top_global, &ems, &top, &f, &f_global,
1671 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1672 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
1673 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1674 * so we free some memory again.
1679 xx = (real *)state->x;
1683 end = mdatoms->homenr;
1685 /* Print to log file */
1686 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1688 do_log = do_ene = do_x = do_f = TRUE;
1690 /* Max number of steps */
1691 number_steps = inputrec->nsteps;
1693 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1695 for (i = start; i < end; i++)
1697 if (mdatoms->cFREEZE)
1699 gf = mdatoms->cFREEZE[i];
1701 for (m = 0; m < DIM; m++)
1703 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1708 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1712 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1717 construct_vsites(vsite, state->x, 1, NULL,
1718 top->idef.iparams, top->idef.il,
1719 fr->ePBC, fr->bMolPBC, cr, state->box);
1722 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1723 /* do_force always puts the charge groups in the box and shifts again
1724 * We do not unshift, so molecules are always whole
1729 evaluate_energy(fplog, cr,
1730 top_global, &ems, top,
1731 inputrec, nrnb, wcycle, gstat,
1732 vsite, constr, fcd, graph, mdatoms, fr,
1733 mu_tot, enerd, vir, pres, -1, TRUE);
1738 /* Copy stuff to the energy bin for easy printing etc. */
1739 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1740 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
1741 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1743 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
1744 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1745 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1749 /* This is the starting energy */
1750 Epot = enerd->term[F_EPOT];
1756 /* Set the initial step.
1757 * since it will be multiplied by the non-normalized search direction
1758 * vector (force vector the first time), we scale it by the
1759 * norm of the force.
1764 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1765 fprintf(stderr, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1766 fprintf(stderr, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1767 fprintf(stderr, "\n");
1768 /* and copy to the log file too... */
1769 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1770 fprintf(fplog, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1771 fprintf(fplog, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1772 fprintf(fplog, "\n");
1776 for (i = 0; i < n; i++)
1780 dx[point][i] = ff[i]; /* Initial search direction */
1788 stepsize = 1.0/fnorm;
1790 /* Start the loop over BFGS steps.
1791 * Each successful step is counted, and we continue until
1792 * we either converge or reach the max number of steps.
1797 /* Set the gradient from the force */
1799 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1802 /* Write coordinates if necessary */
1803 do_x = do_per_step(step, inputrec->nstxout);
1804 do_f = do_per_step(step, inputrec->nstfout);
1809 mdof_flags |= MDOF_X;
1814 mdof_flags |= MDOF_F;
1819 mdof_flags |= MDOF_IMD;
1822 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1823 top_global, step, (real)step, state, state, f, f);
1825 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1827 /* pointer to current direction - point=0 first time here */
1830 /* calculate line gradient */
1831 for (gpa = 0, i = 0; i < n; i++)
1836 /* Calculate minimum allowed stepsize, before the average (norm)
1837 * relative change in coordinate is smaller than precision
1839 for (minstep = 0, i = 0; i < n; i++)
1849 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1851 if (stepsize < minstep)
1857 /* Store old forces and coordinates */
1858 for (i = 0; i < n; i++)
1867 for (i = 0; i < n; i++)
1872 /* Take a step downhill.
1873 * In theory, we should minimize the function along this direction.
1874 * That is quite possible, but it turns out to take 5-10 function evaluations
1875 * for each line. However, we dont really need to find the exact minimum -
1876 * it is much better to start a new BFGS step in a modified direction as soon
1877 * as we are close to it. This will save a lot of energy evaluations.
1879 * In practice, we just try to take a single step.
1880 * If it worked (i.e. lowered the energy), we increase the stepsize but
1881 * the continue straight to the next BFGS step without trying to find any minimum.
1882 * If it didn't work (higher energy), there must be a minimum somewhere between
1883 * the old position and the new one.
1885 * Due to the finite numerical accuracy, it turns out that it is a good idea
1886 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1887 * This leads to lower final energies in the tests I've done. / Erik
1892 c = a + stepsize; /* reference position along line is zero */
1894 /* Check stepsize first. We do not allow displacements
1895 * larger than emstep.
1901 for (i = 0; i < n; i++)
1904 if (delta > maxdelta)
1909 if (maxdelta > inputrec->em_stepsize)
1914 while (maxdelta > inputrec->em_stepsize);
1916 /* Take a trial step */
1917 for (i = 0; i < n; i++)
1919 xc[i] = lastx[i] + c*s[i];
1923 /* Calculate energy for the trial step */
1924 ems.s.x = (rvec *)xc;
1926 evaluate_energy(fplog, cr,
1927 top_global, &ems, top,
1928 inputrec, nrnb, wcycle, gstat,
1929 vsite, constr, fcd, graph, mdatoms, fr,
1930 mu_tot, enerd, vir, pres, step, FALSE);
1933 /* Calc derivative along line */
1934 for (gpc = 0, i = 0; i < n; i++)
1936 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1938 /* Sum the gradient along the line across CPUs */
1941 gmx_sumd(1, &gpc, cr);
1944 /* This is the max amount of increase in energy we tolerate */
1945 tmp = sqrt(GMX_REAL_EPS)*fabs(EpotA);
1947 /* Accept the step if the energy is lower, or if it is not significantly higher
1948 * and the line derivative is still negative.
1950 if (EpotC < EpotA || (gpc < 0 && EpotC < (EpotA+tmp)))
1953 /* Great, we found a better energy. Increase step for next iteration
1954 * if we are still going down, decrease it otherwise
1958 stepsize *= 1.618034; /* The golden section */
1962 stepsize *= 0.618034; /* 1/golden section */
1967 /* New energy is the same or higher. We will have to do some work
1968 * to find a smaller value in the interval. Take smaller step next time!
1971 stepsize *= 0.618034;
1974 /* OK, if we didn't find a lower value we will have to locate one now - there must
1975 * be one in the interval [a=0,c].
1976 * The same thing is valid here, though: Don't spend dozens of iterations to find
1977 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1978 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1980 * I also have a safeguard for potentially really patological functions so we never
1981 * take more than 20 steps before we give up ...
1983 * If we already found a lower value we just skip this step and continue to the update.
1992 /* Select a new trial point.
1993 * If the derivatives at points a & c have different sign we interpolate to zero,
1994 * otherwise just do a bisection.
1997 if (gpa < 0 && gpc > 0)
1999 b = a + gpa*(a-c)/(gpc-gpa);
2006 /* safeguard if interpolation close to machine accuracy causes errors:
2007 * never go outside the interval
2009 if (b <= a || b >= c)
2014 /* Take a trial step */
2015 for (i = 0; i < n; i++)
2017 xb[i] = lastx[i] + b*s[i];
2021 /* Calculate energy for the trial step */
2022 ems.s.x = (rvec *)xb;
2024 evaluate_energy(fplog, cr,
2025 top_global, &ems, top,
2026 inputrec, nrnb, wcycle, gstat,
2027 vsite, constr, fcd, graph, mdatoms, fr,
2028 mu_tot, enerd, vir, pres, step, FALSE);
2033 for (gpb = 0, i = 0; i < n; i++)
2035 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2038 /* Sum the gradient along the line across CPUs */
2041 gmx_sumd(1, &gpb, cr);
2044 /* Keep one of the intervals based on the value of the derivative at the new point */
2047 /* Replace c endpoint with b */
2051 /* swap coord pointers b/c */
2061 /* Replace a endpoint with b */
2065 /* swap coord pointers a/b */
2075 * Stop search as soon as we find a value smaller than the endpoints,
2076 * or if the tolerance is below machine precision.
2077 * Never run more than 20 steps, no matter what.
2081 while ((EpotB > EpotA || EpotB > EpotC) && (nminstep < 20));
2083 if (fabs(EpotB-Epot0) < GMX_REAL_EPS || nminstep >= 20)
2085 /* OK. We couldn't find a significantly lower energy.
2086 * If ncorr==0 this was steepest descent, and then we give up.
2087 * If not, reset memory to restart as steepest descent before quitting.
2099 /* Search in gradient direction */
2100 for (i = 0; i < n; i++)
2102 dx[point][i] = ff[i];
2104 /* Reset stepsize */
2105 stepsize = 1.0/fnorm;
2110 /* Select min energy state of A & C, put the best in xx/ff/Epot
2116 for (i = 0; i < n; i++)
2127 for (i = 0; i < n; i++)
2141 for (i = 0; i < n; i++)
2149 /* Update the memory information, and calculate a new
2150 * approximation of the inverse hessian
2153 /* Have new data in Epot, xx, ff */
2154 if (ncorr < nmaxcorr)
2159 for (i = 0; i < n; i++)
2161 dg[point][i] = lastf[i]-ff[i];
2162 dx[point][i] *= stepsize;
2167 for (i = 0; i < n; i++)
2169 dgdg += dg[point][i]*dg[point][i];
2170 dgdx += dg[point][i]*dx[point][i];
2175 rho[point] = 1.0/dgdx;
2178 if (point >= nmaxcorr)
2184 for (i = 0; i < n; i++)
2191 /* Recursive update. First go back over the memory points */
2192 for (k = 0; k < ncorr; k++)
2201 for (i = 0; i < n; i++)
2203 sq += dx[cp][i]*p[i];
2206 alpha[cp] = rho[cp]*sq;
2208 for (i = 0; i < n; i++)
2210 p[i] -= alpha[cp]*dg[cp][i];
2214 for (i = 0; i < n; i++)
2219 /* And then go forward again */
2220 for (k = 0; k < ncorr; k++)
2223 for (i = 0; i < n; i++)
2225 yr += p[i]*dg[cp][i];
2229 beta = alpha[cp]-beta;
2231 for (i = 0; i < n; i++)
2233 p[i] += beta*dx[cp][i];
2243 for (i = 0; i < n; i++)
2247 dx[point][i] = p[i];
2257 /* Test whether the convergence criterion is met */
2258 get_f_norm_max(cr, &(inputrec->opts), mdatoms, f, &fnorm, &fmax, &nfmax);
2260 /* Print it if necessary */
2265 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2266 step, Epot, fnorm/sqrt(state->natoms), fmax, nfmax+1);
2268 /* Store the new (lower) energies */
2269 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2270 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
2271 NULL, NULL, vir, pres, NULL, mu_tot, constr);
2272 do_log = do_per_step(step, inputrec->nstlog);
2273 do_ene = do_per_step(step, inputrec->nstenergy);
2276 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2278 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2279 do_log ? fplog : NULL, step, step, eprNORMAL,
2280 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2283 /* Send x and E to IMD client, if bIMD is TRUE. */
2284 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state->box, state->x, inputrec, 0, wcycle) && MASTER(cr))
2286 IMD_send_positions(inputrec->imd);
2289 /* Stop when the maximum force lies below tolerance.
2290 * If we have reached machine precision, converged is already set to true.
2293 converged = converged || (fmax < inputrec->em_tol);
2295 } /* End of the loop */
2297 /* IMD cleanup, if bIMD is TRUE. */
2298 IMD_finalize(inputrec->bIMD, inputrec->imd);
2302 step--; /* we never took that last step in this case */
2305 if (fmax > inputrec->em_tol)
2309 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2310 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2315 /* If we printed energy and/or logfile last step (which was the last step)
2316 * we don't have to do it again, but otherwise print the final values.
2318 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2320 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2322 if (!do_ene || !do_log) /* Write final energy file entries */
2324 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2325 !do_log ? fplog : NULL, step, step, eprNORMAL,
2326 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2329 /* Print some stuff... */
2332 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2336 * For accurate normal mode calculation it is imperative that we
2337 * store the last conformation into the full precision binary trajectory.
2339 * However, we should only do it if we did NOT already write this step
2340 * above (which we did if do_x or do_f was true).
2342 do_x = !do_per_step(step, inputrec->nstxout);
2343 do_f = !do_per_step(step, inputrec->nstfout);
2344 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2345 top_global, inputrec, step,
2350 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2351 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2352 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2353 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2355 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2358 finish_em(cr, outf, walltime_accounting, wcycle);
2360 /* To print the actual number of steps we needed somewhere */
2361 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2364 } /* That's all folks */
2367 double do_steep(FILE *fplog, t_commrec *cr,
2368 int nfile, const t_filenm fnm[],
2369 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2370 int gmx_unused nstglobalcomm,
2371 gmx_vsite_t *vsite, gmx_constr_t constr,
2372 int gmx_unused stepout,
2373 t_inputrec *inputrec,
2374 gmx_mtop_t *top_global, t_fcdata *fcd,
2375 t_state *state_global,
2377 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2378 gmx_edsam_t gmx_unused ed,
2380 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2381 gmx_membed_t gmx_unused membed,
2382 real gmx_unused cpt_period, real gmx_unused max_hours,
2383 const char gmx_unused *deviceOptions,
2385 unsigned long gmx_unused Flags,
2386 gmx_walltime_accounting_t walltime_accounting)
2388 const char *SD = "Steepest Descents";
2389 em_state_t *s_min, *s_try;
2391 gmx_localtop_t *top;
2392 gmx_enerdata_t *enerd;
2394 gmx_global_stat_t gstat;
2396 real stepsize, constepsize;
2400 gmx_bool bDone, bAbort, do_x, do_f;
2405 int steps_accepted = 0;
2409 s_min = init_em_state();
2410 s_try = init_em_state();
2412 /* Init em and store the local state in s_try */
2413 init_em(fplog, SD, cr, inputrec,
2414 state_global, top_global, s_try, &top, &f, &f_global,
2415 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2416 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
2418 /* Print to log file */
2419 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2421 /* Set variables for stepsize (in nm). This is the largest
2422 * step that we are going to make in any direction.
2424 ustep = inputrec->em_stepsize;
2427 /* Max number of steps */
2428 nsteps = inputrec->nsteps;
2432 /* Print to the screen */
2433 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2437 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2440 /**** HERE STARTS THE LOOP ****
2441 * count is the counter for the number of steps
2442 * bDone will be TRUE when the minimization has converged
2443 * bAbort will be TRUE when nsteps steps have been performed or when
2444 * the stepsize becomes smaller than is reasonable for machine precision
2449 while (!bDone && !bAbort)
2451 bAbort = (nsteps >= 0) && (count == nsteps);
2453 /* set new coordinates, except for first step */
2456 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2457 s_min, stepsize, s_min->f, s_try,
2458 constr, top, nrnb, wcycle, count);
2461 evaluate_energy(fplog, cr,
2462 top_global, s_try, top,
2463 inputrec, nrnb, wcycle, gstat,
2464 vsite, constr, fcd, graph, mdatoms, fr,
2465 mu_tot, enerd, vir, pres, count, count == 0);
2469 print_ebin_header(fplog, count, count, s_try->s.lambda[efptFEP]);
2474 s_min->epot = s_try->epot + 1;
2477 /* Print it if necessary */
2482 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2483 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2484 (s_try->epot < s_min->epot) ? '\n' : '\r');
2487 if (s_try->epot < s_min->epot)
2489 /* Store the new (lower) energies */
2490 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2491 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2492 s_try->s.box, NULL, NULL, vir, pres, NULL, mu_tot, constr);
2494 /* Prepare IMD energy record, if bIMD is TRUE. */
2495 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2497 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2498 do_per_step(steps_accepted, inputrec->nstdisreout),
2499 do_per_step(steps_accepted, inputrec->nstorireout),
2500 fplog, count, count, eprNORMAL, TRUE,
2501 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2506 /* Now if the new energy is smaller than the previous...
2507 * or if this is the first step!
2508 * or if we did random steps!
2511 if ( (count == 0) || (s_try->epot < s_min->epot) )
2515 /* Test whether the convergence criterion is met... */
2516 bDone = (s_try->fmax < inputrec->em_tol);
2518 /* Copy the arrays for force, positions and energy */
2519 /* The 'Min' array always holds the coords and forces of the minimal
2521 swap_em_state(s_min, s_try);
2527 /* Write to trn, if necessary */
2528 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2529 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2530 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
2531 top_global, inputrec, count,
2532 s_min, state_global, f_global);
2536 /* If energy is not smaller make the step smaller... */
2539 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2541 /* Reload the old state */
2542 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2543 s_min, top, mdatoms, fr, vsite, constr,
2548 /* Determine new step */
2549 stepsize = ustep/s_min->fmax;
2551 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2553 if (count == nsteps || ustep < 1e-12)
2555 if (count == nsteps || ustep < 1e-6)
2560 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2561 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2566 /* Send IMD energies and positions, if bIMD is TRUE. */
2567 if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
2569 IMD_send_positions(inputrec->imd);
2573 } /* End of the loop */
2575 /* IMD cleanup, if bIMD is TRUE. */
2576 IMD_finalize(inputrec->bIMD, inputrec->imd);
2578 /* Print some data... */
2581 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2583 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2584 top_global, inputrec, count,
2585 s_min, state_global, f_global);
2587 fnormn = s_min->fnorm/sqrt(state_global->natoms);
2591 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2592 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2593 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2594 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2597 finish_em(cr, outf, walltime_accounting, wcycle);
2599 /* To print the actual number of steps we needed somewhere */
2600 inputrec->nsteps = count;
2602 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2605 } /* That's all folks */
2608 double do_nm(FILE *fplog, t_commrec *cr,
2609 int nfile, const t_filenm fnm[],
2610 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2611 int gmx_unused nstglobalcomm,
2612 gmx_vsite_t *vsite, gmx_constr_t constr,
2613 int gmx_unused stepout,
2614 t_inputrec *inputrec,
2615 gmx_mtop_t *top_global, t_fcdata *fcd,
2616 t_state *state_global,
2618 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2619 gmx_edsam_t gmx_unused ed,
2621 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2622 gmx_membed_t gmx_unused membed,
2623 real gmx_unused cpt_period, real gmx_unused max_hours,
2624 const char gmx_unused *deviceOptions,
2626 unsigned long gmx_unused Flags,
2627 gmx_walltime_accounting_t walltime_accounting)
2629 const char *NM = "Normal Mode Analysis";
2631 int natoms, atom, d;
2634 gmx_localtop_t *top;
2635 gmx_enerdata_t *enerd;
2637 gmx_global_stat_t gstat;
2639 real t, t0, lambda, lam0;
2644 gmx_bool bSparse; /* use sparse matrix storage format */
2646 gmx_sparsematrix_t * sparse_matrix = NULL;
2647 real * full_matrix = NULL;
2648 em_state_t * state_work;
2650 /* added with respect to mdrun */
2651 int i, j, k, row, col;
2652 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2658 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2661 state_work = init_em_state();
2663 /* Init em and store the local state in state_minimum */
2664 init_em(fplog, NM, cr, inputrec,
2665 state_global, top_global, state_work, &top,
2667 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2668 nfile, fnm, &outf, NULL, imdport, Flags, wcycle);
2670 natoms = top_global->natoms;
2678 "NOTE: This version of Gromacs has been compiled in single precision,\n"
2679 " which MIGHT not be accurate enough for normal mode analysis.\n"
2680 " Gromacs now uses sparse matrix storage, so the memory requirements\n"
2681 " are fairly modest even if you recompile in double precision.\n\n");
2685 /* Check if we can/should use sparse storage format.
2687 * Sparse format is only useful when the Hessian itself is sparse, which it
2688 * will be when we use a cutoff.
2689 * For small systems (n<1000) it is easier to always use full matrix format, though.
2691 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2693 md_print_info(cr, fplog, "Non-cutoff electrostatics used, forcing full Hessian format.\n");
2696 else if (top_global->natoms < 1000)
2698 md_print_info(cr, fplog, "Small system size (N=%d), using full Hessian format.\n", top_global->natoms);
2703 md_print_info(cr, fplog, "Using compressed symmetric sparse Hessian format.\n");
2709 sz = DIM*top_global->natoms;
2711 fprintf(stderr, "Allocating Hessian memory...\n\n");
2715 sparse_matrix = gmx_sparsematrix_init(sz);
2716 sparse_matrix->compressed_symmetric = TRUE;
2720 snew(full_matrix, sz*sz);
2724 /* Initial values */
2725 t0 = inputrec->init_t;
2726 lam0 = inputrec->fepvals->init_lambda;
2734 /* Write start time and temperature */
2735 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2737 /* fudge nr of steps to nr of atoms */
2738 inputrec->nsteps = natoms*2;
2742 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2743 *(top_global->name), (int)inputrec->nsteps);
2746 nnodes = cr->nnodes;
2748 /* Make evaluate_energy do a single node force calculation */
2750 evaluate_energy(fplog, cr,
2751 top_global, state_work, top,
2752 inputrec, nrnb, wcycle, gstat,
2753 vsite, constr, fcd, graph, mdatoms, fr,
2754 mu_tot, enerd, vir, pres, -1, TRUE);
2755 cr->nnodes = nnodes;
2757 /* if forces are not small, warn user */
2758 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, state_work);
2760 md_print_info(cr, fplog, "Maximum force:%12.5e\n", state_work->fmax);
2761 if (state_work->fmax > 1.0e-3)
2763 md_print_info(cr, fplog,
2764 "The force is probably not small enough to "
2765 "ensure that you are at a minimum.\n"
2766 "Be aware that negative eigenvalues may occur\n"
2767 "when the resulting matrix is diagonalized.\n\n");
2770 /***********************************************************
2772 * Loop over all pairs in matrix
2774 * do_force called twice. Once with positive and
2775 * once with negative displacement
2777 ************************************************************/
2779 /* Steps are divided one by one over the nodes */
2780 for (atom = cr->nodeid; atom < natoms; atom += nnodes)
2783 for (d = 0; d < DIM; d++)
2785 x_min = state_work->s.x[atom][d];
2787 state_work->s.x[atom][d] = x_min - der_range;
2789 /* Make evaluate_energy do a single node force calculation */
2791 evaluate_energy(fplog, cr,
2792 top_global, state_work, top,
2793 inputrec, nrnb, wcycle, gstat,
2794 vsite, constr, fcd, graph, mdatoms, fr,
2795 mu_tot, enerd, vir, pres, atom*2, FALSE);
2797 for (i = 0; i < natoms; i++)
2799 copy_rvec(state_work->f[i], fneg[i]);
2802 state_work->s.x[atom][d] = x_min + der_range;
2804 evaluate_energy(fplog, cr,
2805 top_global, state_work, top,
2806 inputrec, nrnb, wcycle, gstat,
2807 vsite, constr, fcd, graph, mdatoms, fr,
2808 mu_tot, enerd, vir, pres, atom*2+1, FALSE);
2809 cr->nnodes = nnodes;
2811 /* x is restored to original */
2812 state_work->s.x[atom][d] = x_min;
2814 for (j = 0; j < natoms; j++)
2816 for (k = 0; (k < DIM); k++)
2819 -(state_work->f[j][k] - fneg[j][k])/(2*der_range);
2827 #define mpi_type MPI_DOUBLE
2829 #define mpi_type MPI_FLOAT
2831 MPI_Send(dfdx[0], natoms*DIM, mpi_type, MASTERNODE(cr), cr->nodeid,
2832 cr->mpi_comm_mygroup);
2837 for (node = 0; (node < nnodes && atom+node < natoms); node++)
2843 MPI_Recv(dfdx[0], natoms*DIM, mpi_type, node, node,
2844 cr->mpi_comm_mygroup, &stat);
2849 row = (atom + node)*DIM + d;
2851 for (j = 0; j < natoms; j++)
2853 for (k = 0; k < DIM; k++)
2859 if (col >= row && dfdx[j][k] != 0.0)
2861 gmx_sparsematrix_increment_value(sparse_matrix,
2862 row, col, dfdx[j][k]);
2867 full_matrix[row*sz+col] = dfdx[j][k];
2874 if (bVerbose && fplog)
2879 /* write progress */
2880 if (MASTER(cr) && bVerbose)
2882 fprintf(stderr, "\rFinished step %d out of %d",
2883 min(atom+nnodes, natoms), natoms);
2890 fprintf(stderr, "\n\nWriting Hessian...\n");
2891 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2894 finish_em(cr, outf, walltime_accounting, wcycle);
2896 walltime_accounting_set_nsteps_done(walltime_accounting, natoms*2);