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/legacyheaders/network.h"
46 #include "gromacs/legacyheaders/nrnb.h"
47 #include "gromacs/legacyheaders/force.h"
48 #include "gromacs/legacyheaders/macros.h"
49 #include "gromacs/legacyheaders/names.h"
50 #include "gromacs/legacyheaders/txtdump.h"
51 #include "gromacs/legacyheaders/typedefs.h"
52 #include "gromacs/legacyheaders/update.h"
53 #include "gromacs/legacyheaders/constr.h"
54 #include "gromacs/legacyheaders/tgroup.h"
55 #include "gromacs/legacyheaders/mdebin.h"
56 #include "gromacs/legacyheaders/vsite.h"
57 #include "gromacs/legacyheaders/force.h"
58 #include "gromacs/legacyheaders/mdrun.h"
59 #include "gromacs/legacyheaders/md_support.h"
60 #include "gromacs/legacyheaders/sim_util.h"
61 #include "gromacs/legacyheaders/domdec.h"
62 #include "gromacs/legacyheaders/mdatoms.h"
63 #include "gromacs/legacyheaders/ns.h"
64 #include "gromacs/topology/mtop_util.h"
65 #include "gromacs/legacyheaders/pme.h"
66 #include "gromacs/legacyheaders/bonded-threading.h"
67 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
68 #include "gromacs/legacyheaders/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/pbcutil/pbc.h"
79 #include "gromacs/timing/wallcycle.h"
80 #include "gromacs/timing/walltime_accounting.h"
81 #include "gromacs/utility/cstringutil.h"
82 #include "gromacs/utility/fatalerror.h"
83 #include "gromacs/utility/smalloc.h"
94 static em_state_t *init_em_state()
100 /* does this need to be here? Should the array be declared differently (staticaly)in the state definition? */
101 snew(ems->s.lambda, efptNR);
106 static void print_em_start(FILE *fplog,
108 gmx_walltime_accounting_t walltime_accounting,
109 gmx_wallcycle_t wcycle,
112 walltime_accounting_start(walltime_accounting);
113 wallcycle_start(wcycle, ewcRUN);
114 print_start(fplog, cr, walltime_accounting, name);
116 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
117 gmx_wallcycle_t wcycle)
119 wallcycle_stop(wcycle, ewcRUN);
121 walltime_accounting_end(walltime_accounting);
124 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
127 fprintf(out, "%s:\n", minimizer);
128 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
129 fprintf(out, " Number of steps = %12d\n", nsteps);
132 static void warn_step(FILE *fp, real ftol, gmx_bool bLastStep, gmx_bool bConstrain)
138 "\nEnergy minimization reached the maximum number "
139 "of steps before the forces reached the requested "
140 "precision Fmax < %g.\n", ftol);
145 "\nEnergy minimization has stopped, but the forces have "
146 "not converged to the requested precision Fmax < %g (which "
147 "may not be possible for your system). It stopped "
148 "because the algorithm tried to make a new step whose size "
149 "was too small, or there was no change in the energy since "
150 "last step. Either way, we regard the minimization as "
151 "converged to within the available machine precision, "
152 "given your starting configuration and EM parameters.\n%s%s",
154 sizeof(real) < sizeof(double) ?
155 "\nDouble precision normally gives you higher accuracy, but "
156 "this is often not needed for preparing to run molecular "
160 "You might need to increase your constraint accuracy, or turn\n"
161 "off constraints altogether (set constraints = none in mdp file)\n" :
164 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
169 static void print_converged(FILE *fp, const char *alg, real ftol,
170 gmx_int64_t count, gmx_bool bDone, gmx_int64_t nsteps,
171 real epot, real fmax, int nfmax, real fnorm)
173 char buf[STEPSTRSIZE];
177 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
178 alg, ftol, gmx_step_str(count, buf));
180 else if (count < nsteps)
182 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
183 "but did not reach the requested Fmax < %g.\n",
184 alg, gmx_step_str(count, buf), ftol);
188 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
189 alg, ftol, gmx_step_str(count, buf));
193 fprintf(fp, "Potential Energy = %21.14e\n", epot);
194 fprintf(fp, "Maximum force = %21.14e on atom %d\n", fmax, nfmax+1);
195 fprintf(fp, "Norm of force = %21.14e\n", fnorm);
197 fprintf(fp, "Potential Energy = %14.7e\n", epot);
198 fprintf(fp, "Maximum force = %14.7e on atom %d\n", fmax, nfmax+1);
199 fprintf(fp, "Norm of force = %14.7e\n", fnorm);
203 static void get_f_norm_max(t_commrec *cr,
204 t_grpopts *opts, t_mdatoms *mdatoms, rvec *f,
205 real *fnorm, real *fmax, int *a_fmax)
208 real fmax2, fmax2_0, fam;
209 int la_max, a_max, start, end, i, m, gf;
211 /* This routine finds the largest force and returns it.
212 * On parallel machines the global max is taken.
219 end = mdatoms->homenr;
220 if (mdatoms->cFREEZE)
222 for (i = start; i < end; i++)
224 gf = mdatoms->cFREEZE[i];
226 for (m = 0; m < DIM; m++)
228 if (!opts->nFreeze[gf][m])
243 for (i = start; i < end; i++)
255 if (la_max >= 0 && DOMAINDECOMP(cr))
257 a_max = cr->dd->gatindex[la_max];
265 snew(sum, 2*cr->nnodes+1);
266 sum[2*cr->nodeid] = fmax2;
267 sum[2*cr->nodeid+1] = a_max;
268 sum[2*cr->nnodes] = fnorm2;
269 gmx_sumd(2*cr->nnodes+1, sum, cr);
270 fnorm2 = sum[2*cr->nnodes];
271 /* Determine the global maximum */
272 for (i = 0; i < cr->nnodes; i++)
274 if (sum[2*i] > fmax2)
277 a_max = (int)(sum[2*i+1] + 0.5);
285 *fnorm = sqrt(fnorm2);
297 static void get_state_f_norm_max(t_commrec *cr,
298 t_grpopts *opts, t_mdatoms *mdatoms,
301 get_f_norm_max(cr, opts, mdatoms, ems->f, &ems->fnorm, &ems->fmax, &ems->a_fmax);
304 void init_em(FILE *fplog, const char *title,
305 t_commrec *cr, t_inputrec *ir,
306 t_state *state_global, gmx_mtop_t *top_global,
307 em_state_t *ems, gmx_localtop_t **top,
308 rvec **f, rvec **f_global,
309 t_nrnb *nrnb, rvec mu_tot,
310 t_forcerec *fr, gmx_enerdata_t **enerd,
311 t_graph **graph, t_mdatoms *mdatoms, gmx_global_stat_t *gstat,
312 gmx_vsite_t *vsite, gmx_constr_t constr,
313 int nfile, const t_filenm fnm[],
314 gmx_mdoutf_t *outf, t_mdebin **mdebin,
315 int imdport, unsigned long gmx_unused Flags)
322 fprintf(fplog, "Initiating %s\n", title);
325 state_global->ngtc = 0;
327 /* Initialize lambda variables */
328 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, NULL);
332 /* Interactive molecular dynamics */
333 init_IMD(ir, cr, top_global, fplog, 1, state_global->x,
334 nfile, fnm, NULL, imdport, Flags);
336 if (DOMAINDECOMP(cr))
338 *top = dd_init_local_top(top_global);
340 dd_init_local_state(cr->dd, state_global, &ems->s);
344 /* Distribute the charge groups over the nodes from the master node */
345 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
346 state_global, top_global, ir,
347 &ems->s, &ems->f, mdatoms, *top,
348 fr, vsite, NULL, constr,
350 dd_store_state(cr->dd, &ems->s);
354 snew(*f_global, top_global->natoms);
364 snew(*f, top_global->natoms);
366 /* Just copy the state */
367 ems->s = *state_global;
368 snew(ems->s.x, ems->s.nalloc);
369 snew(ems->f, ems->s.nalloc);
370 for (i = 0; i < state_global->natoms; i++)
372 copy_rvec(state_global->x[i], ems->s.x[i]);
374 copy_mat(state_global->box, ems->s.box);
376 *top = gmx_mtop_generate_local_top(top_global, ir);
379 forcerec_set_excl_load(fr, *top);
381 setup_bonded_threading(fr, &(*top)->idef);
383 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
385 *graph = mk_graph(fplog, &((*top)->idef), 0, top_global->natoms, FALSE, FALSE);
392 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
393 update_mdatoms(mdatoms, state_global->lambda[efptFEP]);
397 set_vsite_top(vsite, *top, mdatoms, cr);
403 if (ir->eConstrAlg == econtSHAKE &&
404 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
406 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
407 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
410 if (!DOMAINDECOMP(cr))
412 set_constraints(constr, *top, ir, mdatoms, cr);
415 if (!ir->bContinuation)
417 /* Constrain the starting coordinates */
419 constrain(PAR(cr) ? NULL : fplog, TRUE, TRUE, constr, &(*top)->idef,
420 ir, NULL, cr, -1, 0, 1.0, mdatoms,
421 ems->s.x, ems->s.x, NULL, fr->bMolPBC, ems->s.box,
422 ems->s.lambda[efptFEP], &dvdl_constr,
423 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
429 *gstat = global_stat_init(ir);
432 *outf = init_mdoutf(fplog, nfile, fnm, 0, cr, ir, top_global, NULL);
435 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
440 /* Init bin for energy stuff */
441 *mdebin = init_mdebin(mdoutf_get_fp_ene(*outf), top_global, ir, NULL);
445 calc_shifts(ems->s.box, fr->shift_vec);
448 static void finish_em(t_commrec *cr, gmx_mdoutf_t outf,
449 gmx_walltime_accounting_t walltime_accounting,
450 gmx_wallcycle_t wcycle)
452 if (!(cr->duty & DUTY_PME))
454 /* Tell the PME only node to finish */
455 gmx_pme_send_finish(cr);
460 em_time_end(walltime_accounting, wcycle);
463 static void swap_em_state(em_state_t *ems1, em_state_t *ems2)
472 static void copy_em_coords(em_state_t *ems, t_state *state)
476 for (i = 0; (i < state->natoms); i++)
478 copy_rvec(ems->s.x[i], state->x[i]);
482 static void write_em_traj(FILE *fplog, t_commrec *cr,
484 gmx_bool bX, gmx_bool bF, const char *confout,
485 gmx_mtop_t *top_global,
486 t_inputrec *ir, gmx_int64_t step,
488 t_state *state_global, rvec *f_global)
491 gmx_bool bIMDout = FALSE;
494 /* Shall we do IMD output? */
497 bIMDout = do_per_step(step, IMD_get_step(ir->imd->setup));
500 if ((bX || bF || bIMDout || confout != NULL) && !DOMAINDECOMP(cr))
502 copy_em_coords(state, state_global);
509 mdof_flags |= MDOF_X;
513 mdof_flags |= MDOF_F;
516 /* If we want IMD output, set appropriate MDOF flag */
519 mdof_flags |= MDOF_IMD;
522 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
523 top_global, step, (double)step,
524 &state->s, state_global, state->f, f_global);
526 if (confout != NULL && MASTER(cr))
528 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
530 /* Make molecules whole only for confout writing */
531 do_pbc_mtop(fplog, ir->ePBC, state_global->box, top_global,
535 write_sto_conf_mtop(confout,
536 *top_global->name, top_global,
537 state_global->x, NULL, ir->ePBC, state_global->box);
541 static void do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
543 em_state_t *ems1, real a, rvec *f, em_state_t *ems2,
544 gmx_constr_t constr, gmx_localtop_t *top,
545 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
558 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
560 gmx_incons("state mismatch in do_em_step");
563 s2->flags = s1->flags;
565 if (s2->nalloc != s1->nalloc)
567 s2->nalloc = s1->nalloc;
568 srenew(s2->x, s1->nalloc);
569 srenew(ems2->f, s1->nalloc);
570 if (s2->flags & (1<<estCGP))
572 srenew(s2->cg_p, s1->nalloc);
576 s2->natoms = s1->natoms;
577 copy_mat(s1->box, s2->box);
578 /* Copy free energy state */
579 for (i = 0; i < efptNR; i++)
581 s2->lambda[i] = s1->lambda[i];
583 copy_mat(s1->box, s2->box);
591 #pragma omp parallel num_threads(gmx_omp_nthreads_get(emntUpdate))
596 #pragma omp for schedule(static) nowait
597 for (i = start; i < end; i++)
603 for (m = 0; m < DIM; m++)
605 if (ir->opts.nFreeze[gf][m])
611 x2[i][m] = x1[i][m] + a*f[i][m];
616 if (s2->flags & (1<<estCGP))
618 /* Copy the CG p vector */
621 #pragma omp for schedule(static) nowait
622 for (i = start; i < end; i++)
624 copy_rvec(x1[i], x2[i]);
628 if (DOMAINDECOMP(cr))
630 s2->ddp_count = s1->ddp_count;
631 if (s2->cg_gl_nalloc < s1->cg_gl_nalloc)
634 s2->cg_gl_nalloc = s1->cg_gl_nalloc;
635 srenew(s2->cg_gl, s2->cg_gl_nalloc);
638 s2->ncg_gl = s1->ncg_gl;
639 #pragma omp for schedule(static) nowait
640 for (i = 0; i < s2->ncg_gl; i++)
642 s2->cg_gl[i] = s1->cg_gl[i];
644 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
650 wallcycle_start(wcycle, ewcCONSTR);
652 constrain(NULL, TRUE, TRUE, constr, &top->idef,
653 ir, NULL, cr, count, 0, 1.0, md,
654 s1->x, s2->x, NULL, bMolPBC, s2->box,
655 s2->lambda[efptBONDED], &dvdl_constr,
656 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
657 wallcycle_stop(wcycle, ewcCONSTR);
661 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
662 gmx_mtop_t *top_global, t_inputrec *ir,
663 em_state_t *ems, gmx_localtop_t *top,
664 t_mdatoms *mdatoms, t_forcerec *fr,
665 gmx_vsite_t *vsite, gmx_constr_t constr,
666 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
668 /* Repartition the domain decomposition */
669 wallcycle_start(wcycle, ewcDOMDEC);
670 dd_partition_system(fplog, step, cr, FALSE, 1,
671 NULL, top_global, ir,
673 mdatoms, top, fr, vsite, NULL, constr,
674 nrnb, wcycle, FALSE);
675 dd_store_state(cr->dd, &ems->s);
676 wallcycle_stop(wcycle, ewcDOMDEC);
679 static void evaluate_energy(FILE *fplog, t_commrec *cr,
680 gmx_mtop_t *top_global,
681 em_state_t *ems, gmx_localtop_t *top,
682 t_inputrec *inputrec,
683 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
684 gmx_global_stat_t gstat,
685 gmx_vsite_t *vsite, gmx_constr_t constr,
687 t_graph *graph, t_mdatoms *mdatoms,
688 t_forcerec *fr, rvec mu_tot,
689 gmx_enerdata_t *enerd, tensor vir, tensor pres,
690 gmx_int64_t count, gmx_bool bFirst)
695 tensor force_vir, shake_vir, ekin;
696 real dvdl_constr, prescorr, enercorr, dvdlcorr;
699 /* Set the time to the initial time, the time does not change during EM */
700 t = inputrec->init_t;
703 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
705 /* This is the first state or an old state used before the last ns */
711 if (inputrec->nstlist > 0)
715 else if (inputrec->nstlist == -1)
717 nabnsb = natoms_beyond_ns_buffer(inputrec, fr, &top->cgs, NULL, ems->s.x);
720 gmx_sumi(1, &nabnsb, cr);
728 construct_vsites(vsite, ems->s.x, 1, NULL,
729 top->idef.iparams, top->idef.il,
730 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
733 if (DOMAINDECOMP(cr) && bNS)
735 /* Repartition the domain decomposition */
736 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
737 ems, top, mdatoms, fr, vsite, constr,
741 /* Calc force & energy on new trial position */
742 /* do_force always puts the charge groups in the box and shifts again
743 * We do not unshift, so molecules are always whole in congrad.c
745 do_force(fplog, cr, inputrec,
746 count, nrnb, wcycle, top, &top_global->groups,
747 ems->s.box, ems->s.x, &ems->s.hist,
748 ems->f, force_vir, mdatoms, enerd, fcd,
749 ems->s.lambda, graph, fr, vsite, mu_tot, t, NULL, NULL, TRUE,
750 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
751 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
752 (bNS ? GMX_FORCE_NS | GMX_FORCE_DO_LR : 0));
754 /* Clear the unused shake virial and pressure */
755 clear_mat(shake_vir);
758 /* Communicate stuff when parallel */
759 if (PAR(cr) && inputrec->eI != eiNM)
761 wallcycle_start(wcycle, ewcMoveE);
763 global_stat(fplog, gstat, cr, enerd, force_vir, shake_vir, mu_tot,
764 inputrec, NULL, NULL, NULL, 1, &terminate,
765 top_global, &ems->s, FALSE,
771 wallcycle_stop(wcycle, ewcMoveE);
774 /* Calculate long range corrections to pressure and energy */
775 calc_dispcorr(inputrec, fr, top_global->natoms, ems->s.box, ems->s.lambda[efptVDW],
776 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
777 enerd->term[F_DISPCORR] = enercorr;
778 enerd->term[F_EPOT] += enercorr;
779 enerd->term[F_PRES] += prescorr;
780 enerd->term[F_DVDL] += dvdlcorr;
782 ems->epot = enerd->term[F_EPOT];
786 /* Project out the constraint components of the force */
787 wallcycle_start(wcycle, ewcCONSTR);
789 constrain(NULL, FALSE, FALSE, constr, &top->idef,
790 inputrec, NULL, cr, count, 0, 1.0, mdatoms,
791 ems->s.x, ems->f, ems->f, fr->bMolPBC, ems->s.box,
792 ems->s.lambda[efptBONDED], &dvdl_constr,
793 NULL, &shake_vir, nrnb, econqForceDispl, FALSE, 0, 0);
794 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
795 m_add(force_vir, shake_vir, vir);
796 wallcycle_stop(wcycle, ewcCONSTR);
800 copy_mat(force_vir, vir);
804 enerd->term[F_PRES] =
805 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
807 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
809 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
811 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
815 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
817 em_state_t *s_min, em_state_t *s_b)
821 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
823 unsigned char *grpnrFREEZE;
827 fprintf(debug, "Doing reorder_partsum\n");
833 cgs_gl = dd_charge_groups_global(cr->dd);
834 index = cgs_gl->index;
836 /* Collect fm in a global vector fmg.
837 * This conflicts with the spirit of domain decomposition,
838 * but to fully optimize this a much more complicated algorithm is required.
840 snew(fmg, mtop->natoms);
842 ncg = s_min->s.ncg_gl;
843 cg_gl = s_min->s.cg_gl;
845 for (c = 0; c < ncg; c++)
850 for (a = a0; a < a1; a++)
852 copy_rvec(fm[i], fmg[a]);
856 gmx_sum(mtop->natoms*3, fmg[0], cr);
858 /* Now we will determine the part of the sum for the cgs in state s_b */
860 cg_gl = s_b->s.cg_gl;
864 grpnrFREEZE = mtop->groups.grpnr[egcFREEZE];
865 for (c = 0; c < ncg; c++)
870 for (a = a0; a < a1; a++)
872 if (mdatoms->cFREEZE && grpnrFREEZE)
876 for (m = 0; m < DIM; m++)
878 if (!opts->nFreeze[gf][m])
880 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
892 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
894 em_state_t *s_min, em_state_t *s_b)
900 /* This is just the classical Polak-Ribiere calculation of beta;
901 * it looks a bit complicated since we take freeze groups into account,
902 * and might have to sum it in parallel runs.
905 if (!DOMAINDECOMP(cr) ||
906 (s_min->s.ddp_count == cr->dd->ddp_count &&
907 s_b->s.ddp_count == cr->dd->ddp_count))
913 /* This part of code can be incorrect with DD,
914 * since the atom ordering in s_b and s_min might differ.
916 for (i = 0; i < mdatoms->homenr; i++)
918 if (mdatoms->cFREEZE)
920 gf = mdatoms->cFREEZE[i];
922 for (m = 0; m < DIM; m++)
924 if (!opts->nFreeze[gf][m])
926 sum += (fb[i][m] - fm[i][m])*fb[i][m];
933 /* We need to reorder cgs while summing */
934 sum = reorder_partsum(cr, opts, mdatoms, mtop, s_min, s_b);
938 gmx_sumd(1, &sum, cr);
941 return sum/sqr(s_min->fnorm);
944 double do_cg(FILE *fplog, t_commrec *cr,
945 int nfile, const t_filenm fnm[],
946 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
947 int gmx_unused nstglobalcomm,
948 gmx_vsite_t *vsite, gmx_constr_t constr,
949 int gmx_unused stepout,
950 t_inputrec *inputrec,
951 gmx_mtop_t *top_global, t_fcdata *fcd,
952 t_state *state_global,
954 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
955 gmx_edsam_t gmx_unused ed,
957 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
958 gmx_membed_t gmx_unused membed,
959 real gmx_unused cpt_period, real gmx_unused max_hours,
960 const char gmx_unused *deviceOptions,
962 unsigned long gmx_unused Flags,
963 gmx_walltime_accounting_t walltime_accounting)
965 const char *CG = "Polak-Ribiere Conjugate Gradients";
967 em_state_t *s_min, *s_a, *s_b, *s_c;
969 gmx_enerdata_t *enerd;
971 gmx_global_stat_t gstat;
973 rvec *f_global, *p, *sf, *sfm;
974 double gpa, gpb, gpc, tmp, sum[2], minstep;
977 real a, b, c, beta = 0.0;
981 gmx_bool converged, foundlower;
983 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
985 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
987 int i, m, gf, step, nminstep;
992 s_min = init_em_state();
993 s_a = init_em_state();
994 s_b = init_em_state();
995 s_c = init_em_state();
997 /* Init em and store the local state in s_min */
998 init_em(fplog, CG, cr, inputrec,
999 state_global, top_global, s_min, &top, &f, &f_global,
1000 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1001 nfile, fnm, &outf, &mdebin, imdport, Flags);
1003 /* Print to log file */
1004 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1006 /* Max number of steps */
1007 number_steps = inputrec->nsteps;
1011 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1015 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1018 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1019 /* do_force always puts the charge groups in the box and shifts again
1020 * We do not unshift, so molecules are always whole in congrad.c
1022 evaluate_energy(fplog, cr,
1023 top_global, s_min, top,
1024 inputrec, nrnb, wcycle, gstat,
1025 vsite, constr, fcd, graph, mdatoms, fr,
1026 mu_tot, enerd, vir, pres, -1, TRUE);
1031 /* Copy stuff to the energy bin for easy printing etc. */
1032 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1033 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1034 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1036 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1037 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1038 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1042 /* Estimate/guess the initial stepsize */
1043 stepsize = inputrec->em_stepsize/s_min->fnorm;
1047 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1048 s_min->fmax, s_min->a_fmax+1);
1049 fprintf(stderr, " F-Norm = %12.5e\n",
1050 s_min->fnorm/sqrt(state_global->natoms));
1051 fprintf(stderr, "\n");
1052 /* and copy to the log file too... */
1053 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1054 s_min->fmax, s_min->a_fmax+1);
1055 fprintf(fplog, " F-Norm = %12.5e\n",
1056 s_min->fnorm/sqrt(state_global->natoms));
1057 fprintf(fplog, "\n");
1059 /* Start the loop over CG steps.
1060 * Each successful step is counted, and we continue until
1061 * we either converge or reach the max number of steps.
1064 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1067 /* start taking steps in a new direction
1068 * First time we enter the routine, beta=0, and the direction is
1069 * simply the negative gradient.
1072 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1077 for (i = 0; i < mdatoms->homenr; i++)
1079 if (mdatoms->cFREEZE)
1081 gf = mdatoms->cFREEZE[i];
1083 for (m = 0; m < DIM; m++)
1085 if (!inputrec->opts.nFreeze[gf][m])
1087 p[i][m] = sf[i][m] + beta*p[i][m];
1088 gpa -= p[i][m]*sf[i][m];
1089 /* f is negative gradient, thus the sign */
1098 /* Sum the gradient along the line across CPUs */
1101 gmx_sumd(1, &gpa, cr);
1104 /* Calculate the norm of the search vector */
1105 get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1107 /* Just in case stepsize reaches zero due to numerical precision... */
1110 stepsize = inputrec->em_stepsize/pnorm;
1114 * Double check the value of the derivative in the search direction.
1115 * If it is positive it must be due to the old information in the
1116 * CG formula, so just remove that and start over with beta=0.
1117 * This corresponds to a steepest descent step.
1122 step--; /* Don't count this step since we are restarting */
1123 continue; /* Go back to the beginning of the big for-loop */
1126 /* Calculate minimum allowed stepsize, before the average (norm)
1127 * relative change in coordinate is smaller than precision
1130 for (i = 0; i < mdatoms->homenr; i++)
1132 for (m = 0; m < DIM; m++)
1134 tmp = fabs(s_min->s.x[i][m]);
1143 /* Add up from all CPUs */
1146 gmx_sumd(1, &minstep, cr);
1149 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1151 if (stepsize < minstep)
1157 /* Write coordinates if necessary */
1158 do_x = do_per_step(step, inputrec->nstxout);
1159 do_f = do_per_step(step, inputrec->nstfout);
1161 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1162 top_global, inputrec, step,
1163 s_min, state_global, f_global);
1165 /* Take a step downhill.
1166 * In theory, we should minimize the function along this direction.
1167 * That is quite possible, but it turns out to take 5-10 function evaluations
1168 * for each line. However, we dont really need to find the exact minimum -
1169 * it is much better to start a new CG step in a modified direction as soon
1170 * as we are close to it. This will save a lot of energy evaluations.
1172 * In practice, we just try to take a single step.
1173 * If it worked (i.e. lowered the energy), we increase the stepsize but
1174 * the continue straight to the next CG step without trying to find any minimum.
1175 * If it didn't work (higher energy), there must be a minimum somewhere between
1176 * the old position and the new one.
1178 * Due to the finite numerical accuracy, it turns out that it is a good idea
1179 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1180 * This leads to lower final energies in the tests I've done. / Erik
1182 s_a->epot = s_min->epot;
1184 c = a + stepsize; /* reference position along line is zero */
1186 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1188 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1189 s_min, top, mdatoms, fr, vsite, constr,
1193 /* Take a trial step (new coords in s_c) */
1194 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1195 constr, top, nrnb, wcycle, -1);
1198 /* Calculate energy for the trial step */
1199 evaluate_energy(fplog, cr,
1200 top_global, s_c, top,
1201 inputrec, nrnb, wcycle, gstat,
1202 vsite, constr, fcd, graph, mdatoms, fr,
1203 mu_tot, enerd, vir, pres, -1, FALSE);
1205 /* Calc derivative along line */
1209 for (i = 0; i < mdatoms->homenr; i++)
1211 for (m = 0; m < DIM; m++)
1213 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1216 /* Sum the gradient along the line across CPUs */
1219 gmx_sumd(1, &gpc, cr);
1222 /* This is the max amount of increase in energy we tolerate */
1223 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1225 /* Accept the step if the energy is lower, or if it is not significantly higher
1226 * and the line derivative is still negative.
1228 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1231 /* Great, we found a better energy. Increase step for next iteration
1232 * if we are still going down, decrease it otherwise
1236 stepsize *= 1.618034; /* The golden section */
1240 stepsize *= 0.618034; /* 1/golden section */
1245 /* New energy is the same or higher. We will have to do some work
1246 * to find a smaller value in the interval. Take smaller step next time!
1249 stepsize *= 0.618034;
1255 /* OK, if we didn't find a lower value we will have to locate one now - there must
1256 * be one in the interval [a=0,c].
1257 * The same thing is valid here, though: Don't spend dozens of iterations to find
1258 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1259 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1261 * I also have a safeguard for potentially really patological functions so we never
1262 * take more than 20 steps before we give up ...
1264 * If we already found a lower value we just skip this step and continue to the update.
1272 /* Select a new trial point.
1273 * If the derivatives at points a & c have different sign we interpolate to zero,
1274 * otherwise just do a bisection.
1276 if (gpa < 0 && gpc > 0)
1278 b = a + gpa*(a-c)/(gpc-gpa);
1285 /* safeguard if interpolation close to machine accuracy causes errors:
1286 * never go outside the interval
1288 if (b <= a || b >= c)
1293 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1295 /* Reload the old state */
1296 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1297 s_min, top, mdatoms, fr, vsite, constr,
1301 /* Take a trial step to this new point - new coords in s_b */
1302 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1303 constr, top, nrnb, wcycle, -1);
1306 /* Calculate energy for the trial step */
1307 evaluate_energy(fplog, cr,
1308 top_global, s_b, top,
1309 inputrec, nrnb, wcycle, gstat,
1310 vsite, constr, fcd, graph, mdatoms, fr,
1311 mu_tot, enerd, vir, pres, -1, FALSE);
1313 /* p does not change within a step, but since the domain decomposition
1314 * might change, we have to use cg_p of s_b here.
1319 for (i = 0; i < mdatoms->homenr; i++)
1321 for (m = 0; m < DIM; m++)
1323 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1326 /* Sum the gradient along the line across CPUs */
1329 gmx_sumd(1, &gpb, cr);
1334 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1335 s_a->epot, s_b->epot, s_c->epot, gpb);
1338 epot_repl = s_b->epot;
1340 /* Keep one of the intervals based on the value of the derivative at the new point */
1343 /* Replace c endpoint with b */
1344 swap_em_state(s_b, s_c);
1350 /* Replace a endpoint with b */
1351 swap_em_state(s_b, s_a);
1357 * Stop search as soon as we find a value smaller than the endpoints.
1358 * Never run more than 20 steps, no matter what.
1362 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1365 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1368 /* OK. We couldn't find a significantly lower energy.
1369 * If beta==0 this was steepest descent, and then we give up.
1370 * If not, set beta=0 and restart with steepest descent before quitting.
1380 /* Reset memory before giving up */
1386 /* Select min energy state of A & C, put the best in B.
1388 if (s_c->epot < s_a->epot)
1392 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1393 s_c->epot, s_a->epot);
1395 swap_em_state(s_b, s_c);
1403 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1404 s_a->epot, s_c->epot);
1406 swap_em_state(s_b, s_a);
1416 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1419 swap_em_state(s_b, s_c);
1424 /* new search direction */
1425 /* beta = 0 means forget all memory and restart with steepest descents. */
1426 if (nstcg && ((step % nstcg) == 0))
1432 /* s_min->fnorm cannot be zero, because then we would have converged
1436 /* Polak-Ribiere update.
1437 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1439 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1441 /* Limit beta to prevent oscillations */
1442 if (fabs(beta) > 5.0)
1448 /* update positions */
1449 swap_em_state(s_min, s_b);
1452 /* Print it if necessary */
1457 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1458 step, s_min->epot, s_min->fnorm/sqrt(state_global->natoms),
1459 s_min->fmax, s_min->a_fmax+1);
1461 /* Store the new (lower) energies */
1462 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1463 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1464 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1466 do_log = do_per_step(step, inputrec->nstlog);
1467 do_ene = do_per_step(step, inputrec->nstenergy);
1469 /* Prepare IMD energy record, if bIMD is TRUE. */
1470 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1474 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1476 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1477 do_log ? fplog : NULL, step, step, eprNORMAL,
1478 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1481 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1482 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
1484 IMD_send_positions(inputrec->imd);
1487 /* Stop when the maximum force lies below tolerance.
1488 * If we have reached machine precision, converged is already set to true.
1490 converged = converged || (s_min->fmax < inputrec->em_tol);
1492 } /* End of the loop */
1494 /* IMD cleanup, if bIMD is TRUE. */
1495 IMD_finalize(inputrec->bIMD, inputrec->imd);
1499 step--; /* we never took that last step in this case */
1502 if (s_min->fmax > inputrec->em_tol)
1506 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1507 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1514 /* If we printed energy and/or logfile last step (which was the last step)
1515 * we don't have to do it again, but otherwise print the final values.
1519 /* Write final value to log since we didn't do anything the last step */
1520 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1522 if (!do_ene || !do_log)
1524 /* Write final energy file entries */
1525 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1526 !do_log ? fplog : NULL, step, step, eprNORMAL,
1527 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1531 /* Print some stuff... */
1534 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1538 * For accurate normal mode calculation it is imperative that we
1539 * store the last conformation into the full precision binary trajectory.
1541 * However, we should only do it if we did NOT already write this step
1542 * above (which we did if do_x or do_f was true).
1544 do_x = !do_per_step(step, inputrec->nstxout);
1545 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1547 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1548 top_global, inputrec, step,
1549 s_min, state_global, f_global);
1551 fnormn = s_min->fnorm/sqrt(state_global->natoms);
1555 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1556 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1557 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1558 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1560 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1563 finish_em(cr, outf, walltime_accounting, wcycle);
1565 /* To print the actual number of steps we needed somewhere */
1566 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1569 } /* That's all folks */
1572 double do_lbfgs(FILE *fplog, t_commrec *cr,
1573 int nfile, const t_filenm fnm[],
1574 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
1575 int gmx_unused nstglobalcomm,
1576 gmx_vsite_t *vsite, gmx_constr_t constr,
1577 int gmx_unused stepout,
1578 t_inputrec *inputrec,
1579 gmx_mtop_t *top_global, t_fcdata *fcd,
1582 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1583 gmx_edsam_t gmx_unused ed,
1585 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1586 gmx_membed_t gmx_unused membed,
1587 real gmx_unused cpt_period, real gmx_unused max_hours,
1588 const char gmx_unused *deviceOptions,
1590 unsigned long gmx_unused Flags,
1591 gmx_walltime_accounting_t walltime_accounting)
1593 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1595 gmx_localtop_t *top;
1596 gmx_enerdata_t *enerd;
1598 gmx_global_stat_t gstat;
1601 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1602 double stepsize, gpa, gpb, gpc, tmp, minstep;
1603 real *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1604 real *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1605 real a, b, c, maxdelta, delta;
1606 real diag, Epot0, Epot, EpotA, EpotB, EpotC;
1607 real dgdx, dgdg, sq, yr, beta;
1609 gmx_bool converged, first;
1612 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1614 int start, end, number_steps;
1616 int i, k, m, n, nfmax, gf, step;
1623 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1628 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).");
1631 n = 3*state->natoms;
1632 nmaxcorr = inputrec->nbfgscorr;
1634 /* Allocate memory */
1635 /* Use pointers to real so we dont have to loop over both atoms and
1636 * dimensions all the time...
1637 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1638 * that point to the same memory.
1651 snew(rho, nmaxcorr);
1652 snew(alpha, nmaxcorr);
1655 for (i = 0; i < nmaxcorr; i++)
1661 for (i = 0; i < nmaxcorr; i++)
1670 init_em(fplog, LBFGS, cr, inputrec,
1671 state, top_global, &ems, &top, &f, &f_global,
1672 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1673 nfile, fnm, &outf, &mdebin, imdport, Flags);
1674 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1675 * so we free some memory again.
1680 xx = (real *)state->x;
1684 end = mdatoms->homenr;
1686 /* Print to log file */
1687 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1689 do_log = do_ene = do_x = do_f = TRUE;
1691 /* Max number of steps */
1692 number_steps = inputrec->nsteps;
1694 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1696 for (i = start; i < end; i++)
1698 if (mdatoms->cFREEZE)
1700 gf = mdatoms->cFREEZE[i];
1702 for (m = 0; m < DIM; m++)
1704 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1709 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1713 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1718 construct_vsites(vsite, state->x, 1, NULL,
1719 top->idef.iparams, top->idef.il,
1720 fr->ePBC, fr->bMolPBC, cr, state->box);
1723 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1724 /* do_force always puts the charge groups in the box and shifts again
1725 * We do not unshift, so molecules are always whole
1730 evaluate_energy(fplog, cr,
1731 top_global, &ems, top,
1732 inputrec, nrnb, wcycle, gstat,
1733 vsite, constr, fcd, graph, mdatoms, fr,
1734 mu_tot, enerd, vir, pres, -1, TRUE);
1739 /* Copy stuff to the energy bin for easy printing etc. */
1740 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1741 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
1742 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1744 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
1745 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1746 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1750 /* This is the starting energy */
1751 Epot = enerd->term[F_EPOT];
1757 /* Set the initial step.
1758 * since it will be multiplied by the non-normalized search direction
1759 * vector (force vector the first time), we scale it by the
1760 * norm of the force.
1765 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1766 fprintf(stderr, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1767 fprintf(stderr, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1768 fprintf(stderr, "\n");
1769 /* and copy to the log file too... */
1770 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1771 fprintf(fplog, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1772 fprintf(fplog, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1773 fprintf(fplog, "\n");
1777 for (i = 0; i < n; i++)
1781 dx[point][i] = ff[i]; /* Initial search direction */
1789 stepsize = 1.0/fnorm;
1791 /* Start the loop over BFGS steps.
1792 * Each successful step is counted, and we continue until
1793 * we either converge or reach the max number of steps.
1798 /* Set the gradient from the force */
1800 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1803 /* Write coordinates if necessary */
1804 do_x = do_per_step(step, inputrec->nstxout);
1805 do_f = do_per_step(step, inputrec->nstfout);
1810 mdof_flags |= MDOF_X;
1815 mdof_flags |= MDOF_F;
1820 mdof_flags |= MDOF_IMD;
1823 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1824 top_global, step, (real)step, state, state, f, f);
1826 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1828 /* pointer to current direction - point=0 first time here */
1831 /* calculate line gradient */
1832 for (gpa = 0, i = 0; i < n; i++)
1837 /* Calculate minimum allowed stepsize, before the average (norm)
1838 * relative change in coordinate is smaller than precision
1840 for (minstep = 0, i = 0; i < n; i++)
1850 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1852 if (stepsize < minstep)
1858 /* Store old forces and coordinates */
1859 for (i = 0; i < n; i++)
1868 for (i = 0; i < n; i++)
1873 /* Take a step downhill.
1874 * In theory, we should minimize the function along this direction.
1875 * That is quite possible, but it turns out to take 5-10 function evaluations
1876 * for each line. However, we dont really need to find the exact minimum -
1877 * it is much better to start a new BFGS step in a modified direction as soon
1878 * as we are close to it. This will save a lot of energy evaluations.
1880 * In practice, we just try to take a single step.
1881 * If it worked (i.e. lowered the energy), we increase the stepsize but
1882 * the continue straight to the next BFGS step without trying to find any minimum.
1883 * If it didn't work (higher energy), there must be a minimum somewhere between
1884 * the old position and the new one.
1886 * Due to the finite numerical accuracy, it turns out that it is a good idea
1887 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1888 * This leads to lower final energies in the tests I've done. / Erik
1893 c = a + stepsize; /* reference position along line is zero */
1895 /* Check stepsize first. We do not allow displacements
1896 * larger than emstep.
1902 for (i = 0; i < n; i++)
1905 if (delta > maxdelta)
1910 if (maxdelta > inputrec->em_stepsize)
1915 while (maxdelta > inputrec->em_stepsize);
1917 /* Take a trial step */
1918 for (i = 0; i < n; i++)
1920 xc[i] = lastx[i] + c*s[i];
1924 /* Calculate energy for the trial step */
1925 ems.s.x = (rvec *)xc;
1927 evaluate_energy(fplog, cr,
1928 top_global, &ems, top,
1929 inputrec, nrnb, wcycle, gstat,
1930 vsite, constr, fcd, graph, mdatoms, fr,
1931 mu_tot, enerd, vir, pres, step, FALSE);
1934 /* Calc derivative along line */
1935 for (gpc = 0, i = 0; i < n; i++)
1937 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1939 /* Sum the gradient along the line across CPUs */
1942 gmx_sumd(1, &gpc, cr);
1945 /* This is the max amount of increase in energy we tolerate */
1946 tmp = sqrt(GMX_REAL_EPS)*fabs(EpotA);
1948 /* Accept the step if the energy is lower, or if it is not significantly higher
1949 * and the line derivative is still negative.
1951 if (EpotC < EpotA || (gpc < 0 && EpotC < (EpotA+tmp)))
1954 /* Great, we found a better energy. Increase step for next iteration
1955 * if we are still going down, decrease it otherwise
1959 stepsize *= 1.618034; /* The golden section */
1963 stepsize *= 0.618034; /* 1/golden section */
1968 /* New energy is the same or higher. We will have to do some work
1969 * to find a smaller value in the interval. Take smaller step next time!
1972 stepsize *= 0.618034;
1975 /* OK, if we didn't find a lower value we will have to locate one now - there must
1976 * be one in the interval [a=0,c].
1977 * The same thing is valid here, though: Don't spend dozens of iterations to find
1978 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1979 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1981 * I also have a safeguard for potentially really patological functions so we never
1982 * take more than 20 steps before we give up ...
1984 * If we already found a lower value we just skip this step and continue to the update.
1993 /* Select a new trial point.
1994 * If the derivatives at points a & c have different sign we interpolate to zero,
1995 * otherwise just do a bisection.
1998 if (gpa < 0 && gpc > 0)
2000 b = a + gpa*(a-c)/(gpc-gpa);
2007 /* safeguard if interpolation close to machine accuracy causes errors:
2008 * never go outside the interval
2010 if (b <= a || b >= c)
2015 /* Take a trial step */
2016 for (i = 0; i < n; i++)
2018 xb[i] = lastx[i] + b*s[i];
2022 /* Calculate energy for the trial step */
2023 ems.s.x = (rvec *)xb;
2025 evaluate_energy(fplog, cr,
2026 top_global, &ems, top,
2027 inputrec, nrnb, wcycle, gstat,
2028 vsite, constr, fcd, graph, mdatoms, fr,
2029 mu_tot, enerd, vir, pres, step, FALSE);
2034 for (gpb = 0, i = 0; i < n; i++)
2036 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2039 /* Sum the gradient along the line across CPUs */
2042 gmx_sumd(1, &gpb, cr);
2045 /* Keep one of the intervals based on the value of the derivative at the new point */
2048 /* Replace c endpoint with b */
2052 /* swap coord pointers b/c */
2062 /* Replace a endpoint with b */
2066 /* swap coord pointers a/b */
2076 * Stop search as soon as we find a value smaller than the endpoints,
2077 * or if the tolerance is below machine precision.
2078 * Never run more than 20 steps, no matter what.
2082 while ((EpotB > EpotA || EpotB > EpotC) && (nminstep < 20));
2084 if (fabs(EpotB-Epot0) < GMX_REAL_EPS || nminstep >= 20)
2086 /* OK. We couldn't find a significantly lower energy.
2087 * If ncorr==0 this was steepest descent, and then we give up.
2088 * If not, reset memory to restart as steepest descent before quitting.
2100 /* Search in gradient direction */
2101 for (i = 0; i < n; i++)
2103 dx[point][i] = ff[i];
2105 /* Reset stepsize */
2106 stepsize = 1.0/fnorm;
2111 /* Select min energy state of A & C, put the best in xx/ff/Epot
2117 for (i = 0; i < n; i++)
2128 for (i = 0; i < n; i++)
2142 for (i = 0; i < n; i++)
2150 /* Update the memory information, and calculate a new
2151 * approximation of the inverse hessian
2154 /* Have new data in Epot, xx, ff */
2155 if (ncorr < nmaxcorr)
2160 for (i = 0; i < n; i++)
2162 dg[point][i] = lastf[i]-ff[i];
2163 dx[point][i] *= stepsize;
2168 for (i = 0; i < n; i++)
2170 dgdg += dg[point][i]*dg[point][i];
2171 dgdx += dg[point][i]*dx[point][i];
2176 rho[point] = 1.0/dgdx;
2179 if (point >= nmaxcorr)
2185 for (i = 0; i < n; i++)
2192 /* Recursive update. First go back over the memory points */
2193 for (k = 0; k < ncorr; k++)
2202 for (i = 0; i < n; i++)
2204 sq += dx[cp][i]*p[i];
2207 alpha[cp] = rho[cp]*sq;
2209 for (i = 0; i < n; i++)
2211 p[i] -= alpha[cp]*dg[cp][i];
2215 for (i = 0; i < n; i++)
2220 /* And then go forward again */
2221 for (k = 0; k < ncorr; k++)
2224 for (i = 0; i < n; i++)
2226 yr += p[i]*dg[cp][i];
2230 beta = alpha[cp]-beta;
2232 for (i = 0; i < n; i++)
2234 p[i] += beta*dx[cp][i];
2244 for (i = 0; i < n; i++)
2248 dx[point][i] = p[i];
2258 /* Test whether the convergence criterion is met */
2259 get_f_norm_max(cr, &(inputrec->opts), mdatoms, f, &fnorm, &fmax, &nfmax);
2261 /* Print it if necessary */
2266 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2267 step, Epot, fnorm/sqrt(state->natoms), fmax, nfmax+1);
2269 /* Store the new (lower) energies */
2270 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2271 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
2272 NULL, NULL, vir, pres, NULL, mu_tot, constr);
2273 do_log = do_per_step(step, inputrec->nstlog);
2274 do_ene = do_per_step(step, inputrec->nstenergy);
2277 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2279 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2280 do_log ? fplog : NULL, step, step, eprNORMAL,
2281 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2284 /* Send x and E to IMD client, if bIMD is TRUE. */
2285 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state->box, state->x, inputrec, 0, wcycle) && MASTER(cr))
2287 IMD_send_positions(inputrec->imd);
2290 /* Stop when the maximum force lies below tolerance.
2291 * If we have reached machine precision, converged is already set to true.
2294 converged = converged || (fmax < inputrec->em_tol);
2296 } /* End of the loop */
2298 /* IMD cleanup, if bIMD is TRUE. */
2299 IMD_finalize(inputrec->bIMD, inputrec->imd);
2303 step--; /* we never took that last step in this case */
2306 if (fmax > inputrec->em_tol)
2310 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2311 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2316 /* If we printed energy and/or logfile last step (which was the last step)
2317 * we don't have to do it again, but otherwise print the final values.
2319 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2321 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2323 if (!do_ene || !do_log) /* Write final energy file entries */
2325 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2326 !do_log ? fplog : NULL, step, step, eprNORMAL,
2327 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2330 /* Print some stuff... */
2333 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2337 * For accurate normal mode calculation it is imperative that we
2338 * store the last conformation into the full precision binary trajectory.
2340 * However, we should only do it if we did NOT already write this step
2341 * above (which we did if do_x or do_f was true).
2343 do_x = !do_per_step(step, inputrec->nstxout);
2344 do_f = !do_per_step(step, inputrec->nstfout);
2345 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2346 top_global, inputrec, step,
2351 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2352 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2353 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2354 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2356 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2359 finish_em(cr, outf, walltime_accounting, wcycle);
2361 /* To print the actual number of steps we needed somewhere */
2362 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2365 } /* That's all folks */
2368 double do_steep(FILE *fplog, t_commrec *cr,
2369 int nfile, const t_filenm fnm[],
2370 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2371 int gmx_unused nstglobalcomm,
2372 gmx_vsite_t *vsite, gmx_constr_t constr,
2373 int gmx_unused stepout,
2374 t_inputrec *inputrec,
2375 gmx_mtop_t *top_global, t_fcdata *fcd,
2376 t_state *state_global,
2378 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2379 gmx_edsam_t gmx_unused ed,
2381 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2382 gmx_membed_t gmx_unused membed,
2383 real gmx_unused cpt_period, real gmx_unused max_hours,
2384 const char gmx_unused *deviceOptions,
2386 unsigned long gmx_unused Flags,
2387 gmx_walltime_accounting_t walltime_accounting)
2389 const char *SD = "Steepest Descents";
2390 em_state_t *s_min, *s_try;
2392 gmx_localtop_t *top;
2393 gmx_enerdata_t *enerd;
2395 gmx_global_stat_t gstat;
2397 real stepsize, constepsize;
2401 gmx_bool bDone, bAbort, do_x, do_f;
2406 int steps_accepted = 0;
2410 s_min = init_em_state();
2411 s_try = init_em_state();
2413 /* Init em and store the local state in s_try */
2414 init_em(fplog, SD, cr, inputrec,
2415 state_global, top_global, s_try, &top, &f, &f_global,
2416 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2417 nfile, fnm, &outf, &mdebin, imdport, Flags);
2419 /* Print to log file */
2420 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2422 /* Set variables for stepsize (in nm). This is the largest
2423 * step that we are going to make in any direction.
2425 ustep = inputrec->em_stepsize;
2428 /* Max number of steps */
2429 nsteps = inputrec->nsteps;
2433 /* Print to the screen */
2434 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2438 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2441 /**** HERE STARTS THE LOOP ****
2442 * count is the counter for the number of steps
2443 * bDone will be TRUE when the minimization has converged
2444 * bAbort will be TRUE when nsteps steps have been performed or when
2445 * the stepsize becomes smaller than is reasonable for machine precision
2450 while (!bDone && !bAbort)
2452 bAbort = (nsteps >= 0) && (count == nsteps);
2454 /* set new coordinates, except for first step */
2457 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2458 s_min, stepsize, s_min->f, s_try,
2459 constr, top, nrnb, wcycle, count);
2462 evaluate_energy(fplog, cr,
2463 top_global, s_try, top,
2464 inputrec, nrnb, wcycle, gstat,
2465 vsite, constr, fcd, graph, mdatoms, fr,
2466 mu_tot, enerd, vir, pres, count, count == 0);
2470 print_ebin_header(fplog, count, count, s_try->s.lambda[efptFEP]);
2475 s_min->epot = s_try->epot + 1;
2478 /* Print it if necessary */
2483 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2484 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2485 (s_try->epot < s_min->epot) ? '\n' : '\r');
2488 if (s_try->epot < s_min->epot)
2490 /* Store the new (lower) energies */
2491 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2492 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2493 s_try->s.box, NULL, NULL, vir, pres, NULL, mu_tot, constr);
2495 /* Prepare IMD energy record, if bIMD is TRUE. */
2496 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2498 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2499 do_per_step(steps_accepted, inputrec->nstdisreout),
2500 do_per_step(steps_accepted, inputrec->nstorireout),
2501 fplog, count, count, eprNORMAL, TRUE,
2502 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2507 /* Now if the new energy is smaller than the previous...
2508 * or if this is the first step!
2509 * or if we did random steps!
2512 if ( (count == 0) || (s_try->epot < s_min->epot) )
2516 /* Test whether the convergence criterion is met... */
2517 bDone = (s_try->fmax < inputrec->em_tol);
2519 /* Copy the arrays for force, positions and energy */
2520 /* The 'Min' array always holds the coords and forces of the minimal
2522 swap_em_state(s_min, s_try);
2528 /* Write to trn, if necessary */
2529 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2530 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2531 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
2532 top_global, inputrec, count,
2533 s_min, state_global, f_global);
2537 /* If energy is not smaller make the step smaller... */
2540 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2542 /* Reload the old state */
2543 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2544 s_min, top, mdatoms, fr, vsite, constr,
2549 /* Determine new step */
2550 stepsize = ustep/s_min->fmax;
2552 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2554 if (count == nsteps || ustep < 1e-12)
2556 if (count == nsteps || ustep < 1e-6)
2561 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2562 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2567 /* Send IMD energies and positions, if bIMD is TRUE. */
2568 if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
2570 IMD_send_positions(inputrec->imd);
2574 } /* End of the loop */
2576 /* IMD cleanup, if bIMD is TRUE. */
2577 IMD_finalize(inputrec->bIMD, inputrec->imd);
2579 /* Print some data... */
2582 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2584 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2585 top_global, inputrec, count,
2586 s_min, state_global, f_global);
2588 fnormn = s_min->fnorm/sqrt(state_global->natoms);
2592 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2593 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2594 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2595 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2598 finish_em(cr, outf, walltime_accounting, wcycle);
2600 /* To print the actual number of steps we needed somewhere */
2601 inputrec->nsteps = count;
2603 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2606 } /* That's all folks */
2609 double do_nm(FILE *fplog, t_commrec *cr,
2610 int nfile, const t_filenm fnm[],
2611 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2612 int gmx_unused nstglobalcomm,
2613 gmx_vsite_t *vsite, gmx_constr_t constr,
2614 int gmx_unused stepout,
2615 t_inputrec *inputrec,
2616 gmx_mtop_t *top_global, t_fcdata *fcd,
2617 t_state *state_global,
2619 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2620 gmx_edsam_t gmx_unused ed,
2622 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2623 gmx_membed_t gmx_unused membed,
2624 real gmx_unused cpt_period, real gmx_unused max_hours,
2625 const char gmx_unused *deviceOptions,
2627 unsigned long gmx_unused Flags,
2628 gmx_walltime_accounting_t walltime_accounting)
2630 const char *NM = "Normal Mode Analysis";
2632 int natoms, atom, d;
2635 gmx_localtop_t *top;
2636 gmx_enerdata_t *enerd;
2638 gmx_global_stat_t gstat;
2640 real t, t0, lambda, lam0;
2645 gmx_bool bSparse; /* use sparse matrix storage format */
2647 gmx_sparsematrix_t * sparse_matrix = NULL;
2648 real * full_matrix = NULL;
2649 em_state_t * state_work;
2651 /* added with respect to mdrun */
2652 int i, j, k, row, col;
2653 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2659 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2662 state_work = init_em_state();
2664 /* Init em and store the local state in state_minimum */
2665 init_em(fplog, NM, cr, inputrec,
2666 state_global, top_global, state_work, &top,
2668 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2669 nfile, fnm, &outf, NULL, imdport, Flags);
2671 natoms = top_global->natoms;
2679 "NOTE: This version of Gromacs has been compiled in single precision,\n"
2680 " which MIGHT not be accurate enough for normal mode analysis.\n"
2681 " Gromacs now uses sparse matrix storage, so the memory requirements\n"
2682 " are fairly modest even if you recompile in double precision.\n\n");
2686 /* Check if we can/should use sparse storage format.
2688 * Sparse format is only useful when the Hessian itself is sparse, which it
2689 * will be when we use a cutoff.
2690 * For small systems (n<1000) it is easier to always use full matrix format, though.
2692 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2694 md_print_info(cr, fplog, "Non-cutoff electrostatics used, forcing full Hessian format.\n");
2697 else if (top_global->natoms < 1000)
2699 md_print_info(cr, fplog, "Small system size (N=%d), using full Hessian format.\n", top_global->natoms);
2704 md_print_info(cr, fplog, "Using compressed symmetric sparse Hessian format.\n");
2710 sz = DIM*top_global->natoms;
2712 fprintf(stderr, "Allocating Hessian memory...\n\n");
2716 sparse_matrix = gmx_sparsematrix_init(sz);
2717 sparse_matrix->compressed_symmetric = TRUE;
2721 snew(full_matrix, sz*sz);
2725 /* Initial values */
2726 t0 = inputrec->init_t;
2727 lam0 = inputrec->fepvals->init_lambda;
2735 /* Write start time and temperature */
2736 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2738 /* fudge nr of steps to nr of atoms */
2739 inputrec->nsteps = natoms*2;
2743 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2744 *(top_global->name), (int)inputrec->nsteps);
2747 nnodes = cr->nnodes;
2749 /* Make evaluate_energy do a single node force calculation */
2751 evaluate_energy(fplog, cr,
2752 top_global, state_work, top,
2753 inputrec, nrnb, wcycle, gstat,
2754 vsite, constr, fcd, graph, mdatoms, fr,
2755 mu_tot, enerd, vir, pres, -1, TRUE);
2756 cr->nnodes = nnodes;
2758 /* if forces are not small, warn user */
2759 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, state_work);
2761 md_print_info(cr, fplog, "Maximum force:%12.5e\n", state_work->fmax);
2762 if (state_work->fmax > 1.0e-3)
2764 md_print_info(cr, fplog,
2765 "The force is probably not small enough to "
2766 "ensure that you are at a minimum.\n"
2767 "Be aware that negative eigenvalues may occur\n"
2768 "when the resulting matrix is diagonalized.\n\n");
2771 /***********************************************************
2773 * Loop over all pairs in matrix
2775 * do_force called twice. Once with positive and
2776 * once with negative displacement
2778 ************************************************************/
2780 /* Steps are divided one by one over the nodes */
2781 for (atom = cr->nodeid; atom < natoms; atom += nnodes)
2784 for (d = 0; d < DIM; d++)
2786 x_min = state_work->s.x[atom][d];
2788 state_work->s.x[atom][d] = x_min - der_range;
2790 /* Make evaluate_energy do a single node force calculation */
2792 evaluate_energy(fplog, cr,
2793 top_global, state_work, top,
2794 inputrec, nrnb, wcycle, gstat,
2795 vsite, constr, fcd, graph, mdatoms, fr,
2796 mu_tot, enerd, vir, pres, atom*2, FALSE);
2798 for (i = 0; i < natoms; i++)
2800 copy_rvec(state_work->f[i], fneg[i]);
2803 state_work->s.x[atom][d] = x_min + der_range;
2805 evaluate_energy(fplog, cr,
2806 top_global, state_work, top,
2807 inputrec, nrnb, wcycle, gstat,
2808 vsite, constr, fcd, graph, mdatoms, fr,
2809 mu_tot, enerd, vir, pres, atom*2+1, FALSE);
2810 cr->nnodes = nnodes;
2812 /* x is restored to original */
2813 state_work->s.x[atom][d] = x_min;
2815 for (j = 0; j < natoms; j++)
2817 for (k = 0; (k < DIM); k++)
2820 -(state_work->f[j][k] - fneg[j][k])/(2*der_range);
2828 #define mpi_type MPI_DOUBLE
2830 #define mpi_type MPI_FLOAT
2832 MPI_Send(dfdx[0], natoms*DIM, mpi_type, MASTERNODE(cr), cr->nodeid,
2833 cr->mpi_comm_mygroup);
2838 for (node = 0; (node < nnodes && atom+node < natoms); node++)
2844 MPI_Recv(dfdx[0], natoms*DIM, mpi_type, node, node,
2845 cr->mpi_comm_mygroup, &stat);
2850 row = (atom + node)*DIM + d;
2852 for (j = 0; j < natoms; j++)
2854 for (k = 0; k < DIM; k++)
2860 if (col >= row && dfdx[j][k] != 0.0)
2862 gmx_sparsematrix_increment_value(sparse_matrix,
2863 row, col, dfdx[j][k]);
2868 full_matrix[row*sz+col] = dfdx[j][k];
2875 if (bVerbose && fplog)
2880 /* write progress */
2881 if (MASTER(cr) && bVerbose)
2883 fprintf(stderr, "\rFinished step %d out of %d",
2884 min(atom+nnodes, natoms), natoms);
2891 fprintf(stderr, "\n\nWriting Hessian...\n");
2892 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2895 finish_em(cr, outf, walltime_accounting, wcycle);
2897 walltime_accounting_set_nsteps_done(walltime_accounting, natoms*2);