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,
553 int nthreads gmx_unused;
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 nthreads = gmx_omp_nthreads_get(emntUpdate);
592 #pragma omp parallel num_threads(nthreads)
597 #pragma omp for schedule(static) nowait
598 for (i = start; i < end; i++)
604 for (m = 0; m < DIM; m++)
606 if (ir->opts.nFreeze[gf][m])
612 x2[i][m] = x1[i][m] + a*f[i][m];
617 if (s2->flags & (1<<estCGP))
619 /* Copy the CG p vector */
622 #pragma omp for schedule(static) nowait
623 for (i = start; i < end; i++)
625 copy_rvec(x1[i], x2[i]);
629 if (DOMAINDECOMP(cr))
631 s2->ddp_count = s1->ddp_count;
632 if (s2->cg_gl_nalloc < s1->cg_gl_nalloc)
635 s2->cg_gl_nalloc = s1->cg_gl_nalloc;
636 srenew(s2->cg_gl, s2->cg_gl_nalloc);
639 s2->ncg_gl = s1->ncg_gl;
640 #pragma omp for schedule(static) nowait
641 for (i = 0; i < s2->ncg_gl; i++)
643 s2->cg_gl[i] = s1->cg_gl[i];
645 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
651 wallcycle_start(wcycle, ewcCONSTR);
653 constrain(NULL, TRUE, TRUE, constr, &top->idef,
654 ir, NULL, cr, count, 0, 1.0, md,
655 s1->x, s2->x, NULL, bMolPBC, s2->box,
656 s2->lambda[efptBONDED], &dvdl_constr,
657 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
658 wallcycle_stop(wcycle, ewcCONSTR);
662 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
663 gmx_mtop_t *top_global, t_inputrec *ir,
664 em_state_t *ems, gmx_localtop_t *top,
665 t_mdatoms *mdatoms, t_forcerec *fr,
666 gmx_vsite_t *vsite, gmx_constr_t constr,
667 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
669 /* Repartition the domain decomposition */
670 wallcycle_start(wcycle, ewcDOMDEC);
671 dd_partition_system(fplog, step, cr, FALSE, 1,
672 NULL, top_global, ir,
674 mdatoms, top, fr, vsite, NULL, constr,
675 nrnb, wcycle, FALSE);
676 dd_store_state(cr->dd, &ems->s);
677 wallcycle_stop(wcycle, ewcDOMDEC);
680 static void evaluate_energy(FILE *fplog, t_commrec *cr,
681 gmx_mtop_t *top_global,
682 em_state_t *ems, gmx_localtop_t *top,
683 t_inputrec *inputrec,
684 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
685 gmx_global_stat_t gstat,
686 gmx_vsite_t *vsite, gmx_constr_t constr,
688 t_graph *graph, t_mdatoms *mdatoms,
689 t_forcerec *fr, rvec mu_tot,
690 gmx_enerdata_t *enerd, tensor vir, tensor pres,
691 gmx_int64_t count, gmx_bool bFirst)
696 tensor force_vir, shake_vir, ekin;
697 real dvdl_constr, prescorr, enercorr, dvdlcorr;
700 /* Set the time to the initial time, the time does not change during EM */
701 t = inputrec->init_t;
704 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
706 /* This is the first state or an old state used before the last ns */
712 if (inputrec->nstlist > 0)
716 else if (inputrec->nstlist == -1)
718 nabnsb = natoms_beyond_ns_buffer(inputrec, fr, &top->cgs, NULL, ems->s.x);
721 gmx_sumi(1, &nabnsb, cr);
729 construct_vsites(vsite, ems->s.x, 1, NULL,
730 top->idef.iparams, top->idef.il,
731 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
734 if (DOMAINDECOMP(cr) && bNS)
736 /* Repartition the domain decomposition */
737 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
738 ems, top, mdatoms, fr, vsite, constr,
742 /* Calc force & energy on new trial position */
743 /* do_force always puts the charge groups in the box and shifts again
744 * We do not unshift, so molecules are always whole in congrad.c
746 do_force(fplog, cr, inputrec,
747 count, nrnb, wcycle, top, &top_global->groups,
748 ems->s.box, ems->s.x, &ems->s.hist,
749 ems->f, force_vir, mdatoms, enerd, fcd,
750 ems->s.lambda, graph, fr, vsite, mu_tot, t, NULL, NULL, TRUE,
751 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
752 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
753 (bNS ? GMX_FORCE_NS | GMX_FORCE_DO_LR : 0));
755 /* Clear the unused shake virial and pressure */
756 clear_mat(shake_vir);
759 /* Communicate stuff when parallel */
760 if (PAR(cr) && inputrec->eI != eiNM)
762 wallcycle_start(wcycle, ewcMoveE);
764 global_stat(fplog, gstat, cr, enerd, force_vir, shake_vir, mu_tot,
765 inputrec, NULL, NULL, NULL, 1, &terminate,
766 top_global, &ems->s, FALSE,
772 wallcycle_stop(wcycle, ewcMoveE);
775 /* Calculate long range corrections to pressure and energy */
776 calc_dispcorr(inputrec, fr, top_global->natoms, ems->s.box, ems->s.lambda[efptVDW],
777 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
778 enerd->term[F_DISPCORR] = enercorr;
779 enerd->term[F_EPOT] += enercorr;
780 enerd->term[F_PRES] += prescorr;
781 enerd->term[F_DVDL] += dvdlcorr;
783 ems->epot = enerd->term[F_EPOT];
787 /* Project out the constraint components of the force */
788 wallcycle_start(wcycle, ewcCONSTR);
790 constrain(NULL, FALSE, FALSE, constr, &top->idef,
791 inputrec, NULL, cr, count, 0, 1.0, mdatoms,
792 ems->s.x, ems->f, ems->f, fr->bMolPBC, ems->s.box,
793 ems->s.lambda[efptBONDED], &dvdl_constr,
794 NULL, &shake_vir, nrnb, econqForceDispl, FALSE, 0, 0);
795 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
796 m_add(force_vir, shake_vir, vir);
797 wallcycle_stop(wcycle, ewcCONSTR);
801 copy_mat(force_vir, vir);
805 enerd->term[F_PRES] =
806 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
808 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
810 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
812 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
816 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
818 em_state_t *s_min, em_state_t *s_b)
822 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
824 unsigned char *grpnrFREEZE;
828 fprintf(debug, "Doing reorder_partsum\n");
834 cgs_gl = dd_charge_groups_global(cr->dd);
835 index = cgs_gl->index;
837 /* Collect fm in a global vector fmg.
838 * This conflicts with the spirit of domain decomposition,
839 * but to fully optimize this a much more complicated algorithm is required.
841 snew(fmg, mtop->natoms);
843 ncg = s_min->s.ncg_gl;
844 cg_gl = s_min->s.cg_gl;
846 for (c = 0; c < ncg; c++)
851 for (a = a0; a < a1; a++)
853 copy_rvec(fm[i], fmg[a]);
857 gmx_sum(mtop->natoms*3, fmg[0], cr);
859 /* Now we will determine the part of the sum for the cgs in state s_b */
861 cg_gl = s_b->s.cg_gl;
865 grpnrFREEZE = mtop->groups.grpnr[egcFREEZE];
866 for (c = 0; c < ncg; c++)
871 for (a = a0; a < a1; a++)
873 if (mdatoms->cFREEZE && grpnrFREEZE)
877 for (m = 0; m < DIM; m++)
879 if (!opts->nFreeze[gf][m])
881 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
893 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
895 em_state_t *s_min, em_state_t *s_b)
901 /* This is just the classical Polak-Ribiere calculation of beta;
902 * it looks a bit complicated since we take freeze groups into account,
903 * and might have to sum it in parallel runs.
906 if (!DOMAINDECOMP(cr) ||
907 (s_min->s.ddp_count == cr->dd->ddp_count &&
908 s_b->s.ddp_count == cr->dd->ddp_count))
914 /* This part of code can be incorrect with DD,
915 * since the atom ordering in s_b and s_min might differ.
917 for (i = 0; i < mdatoms->homenr; i++)
919 if (mdatoms->cFREEZE)
921 gf = mdatoms->cFREEZE[i];
923 for (m = 0; m < DIM; m++)
925 if (!opts->nFreeze[gf][m])
927 sum += (fb[i][m] - fm[i][m])*fb[i][m];
934 /* We need to reorder cgs while summing */
935 sum = reorder_partsum(cr, opts, mdatoms, mtop, s_min, s_b);
939 gmx_sumd(1, &sum, cr);
942 return sum/sqr(s_min->fnorm);
945 double do_cg(FILE *fplog, t_commrec *cr,
946 int nfile, const t_filenm fnm[],
947 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
948 int gmx_unused nstglobalcomm,
949 gmx_vsite_t *vsite, gmx_constr_t constr,
950 int gmx_unused stepout,
951 t_inputrec *inputrec,
952 gmx_mtop_t *top_global, t_fcdata *fcd,
953 t_state *state_global,
955 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
956 gmx_edsam_t gmx_unused ed,
958 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
959 gmx_membed_t gmx_unused membed,
960 real gmx_unused cpt_period, real gmx_unused max_hours,
961 const char gmx_unused *deviceOptions,
963 unsigned long gmx_unused Flags,
964 gmx_walltime_accounting_t walltime_accounting)
966 const char *CG = "Polak-Ribiere Conjugate Gradients";
968 em_state_t *s_min, *s_a, *s_b, *s_c;
970 gmx_enerdata_t *enerd;
972 gmx_global_stat_t gstat;
974 rvec *f_global, *p, *sf, *sfm;
975 double gpa, gpb, gpc, tmp, sum[2], minstep;
978 real a, b, c, beta = 0.0;
982 gmx_bool converged, foundlower;
984 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
986 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
988 int i, m, gf, step, nminstep;
993 s_min = init_em_state();
994 s_a = init_em_state();
995 s_b = init_em_state();
996 s_c = init_em_state();
998 /* Init em and store the local state in s_min */
999 init_em(fplog, CG, cr, inputrec,
1000 state_global, top_global, s_min, &top, &f, &f_global,
1001 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1002 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
1004 /* Print to log file */
1005 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1007 /* Max number of steps */
1008 number_steps = inputrec->nsteps;
1012 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1016 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1019 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1020 /* do_force always puts the charge groups in the box and shifts again
1021 * We do not unshift, so molecules are always whole in congrad.c
1023 evaluate_energy(fplog, cr,
1024 top_global, s_min, top,
1025 inputrec, nrnb, wcycle, gstat,
1026 vsite, constr, fcd, graph, mdatoms, fr,
1027 mu_tot, enerd, vir, pres, -1, TRUE);
1032 /* Copy stuff to the energy bin for easy printing etc. */
1033 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1034 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1035 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1037 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1038 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1039 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1043 /* Estimate/guess the initial stepsize */
1044 stepsize = inputrec->em_stepsize/s_min->fnorm;
1048 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1049 s_min->fmax, s_min->a_fmax+1);
1050 fprintf(stderr, " F-Norm = %12.5e\n",
1051 s_min->fnorm/sqrt(state_global->natoms));
1052 fprintf(stderr, "\n");
1053 /* and copy to the log file too... */
1054 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1055 s_min->fmax, s_min->a_fmax+1);
1056 fprintf(fplog, " F-Norm = %12.5e\n",
1057 s_min->fnorm/sqrt(state_global->natoms));
1058 fprintf(fplog, "\n");
1060 /* Start the loop over CG steps.
1061 * Each successful step is counted, and we continue until
1062 * we either converge or reach the max number of steps.
1065 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1068 /* start taking steps in a new direction
1069 * First time we enter the routine, beta=0, and the direction is
1070 * simply the negative gradient.
1073 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1078 for (i = 0; i < mdatoms->homenr; i++)
1080 if (mdatoms->cFREEZE)
1082 gf = mdatoms->cFREEZE[i];
1084 for (m = 0; m < DIM; m++)
1086 if (!inputrec->opts.nFreeze[gf][m])
1088 p[i][m] = sf[i][m] + beta*p[i][m];
1089 gpa -= p[i][m]*sf[i][m];
1090 /* f is negative gradient, thus the sign */
1099 /* Sum the gradient along the line across CPUs */
1102 gmx_sumd(1, &gpa, cr);
1105 /* Calculate the norm of the search vector */
1106 get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1108 /* Just in case stepsize reaches zero due to numerical precision... */
1111 stepsize = inputrec->em_stepsize/pnorm;
1115 * Double check the value of the derivative in the search direction.
1116 * If it is positive it must be due to the old information in the
1117 * CG formula, so just remove that and start over with beta=0.
1118 * This corresponds to a steepest descent step.
1123 step--; /* Don't count this step since we are restarting */
1124 continue; /* Go back to the beginning of the big for-loop */
1127 /* Calculate minimum allowed stepsize, before the average (norm)
1128 * relative change in coordinate is smaller than precision
1131 for (i = 0; i < mdatoms->homenr; i++)
1133 for (m = 0; m < DIM; m++)
1135 tmp = fabs(s_min->s.x[i][m]);
1144 /* Add up from all CPUs */
1147 gmx_sumd(1, &minstep, cr);
1150 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1152 if (stepsize < minstep)
1158 /* Write coordinates if necessary */
1159 do_x = do_per_step(step, inputrec->nstxout);
1160 do_f = do_per_step(step, inputrec->nstfout);
1162 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1163 top_global, inputrec, step,
1164 s_min, state_global, f_global);
1166 /* Take a step downhill.
1167 * In theory, we should minimize the function along this direction.
1168 * That is quite possible, but it turns out to take 5-10 function evaluations
1169 * for each line. However, we dont really need to find the exact minimum -
1170 * it is much better to start a new CG step in a modified direction as soon
1171 * as we are close to it. This will save a lot of energy evaluations.
1173 * In practice, we just try to take a single step.
1174 * If it worked (i.e. lowered the energy), we increase the stepsize but
1175 * the continue straight to the next CG step without trying to find any minimum.
1176 * If it didn't work (higher energy), there must be a minimum somewhere between
1177 * the old position and the new one.
1179 * Due to the finite numerical accuracy, it turns out that it is a good idea
1180 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1181 * This leads to lower final energies in the tests I've done. / Erik
1183 s_a->epot = s_min->epot;
1185 c = a + stepsize; /* reference position along line is zero */
1187 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1189 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1190 s_min, top, mdatoms, fr, vsite, constr,
1194 /* Take a trial step (new coords in s_c) */
1195 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1196 constr, top, nrnb, wcycle, -1);
1199 /* Calculate energy for the trial step */
1200 evaluate_energy(fplog, cr,
1201 top_global, s_c, top,
1202 inputrec, nrnb, wcycle, gstat,
1203 vsite, constr, fcd, graph, mdatoms, fr,
1204 mu_tot, enerd, vir, pres, -1, FALSE);
1206 /* Calc derivative along line */
1210 for (i = 0; i < mdatoms->homenr; i++)
1212 for (m = 0; m < DIM; m++)
1214 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1217 /* Sum the gradient along the line across CPUs */
1220 gmx_sumd(1, &gpc, cr);
1223 /* This is the max amount of increase in energy we tolerate */
1224 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1226 /* Accept the step if the energy is lower, or if it is not significantly higher
1227 * and the line derivative is still negative.
1229 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1232 /* Great, we found a better energy. Increase step for next iteration
1233 * if we are still going down, decrease it otherwise
1237 stepsize *= 1.618034; /* The golden section */
1241 stepsize *= 0.618034; /* 1/golden section */
1246 /* New energy is the same or higher. We will have to do some work
1247 * to find a smaller value in the interval. Take smaller step next time!
1250 stepsize *= 0.618034;
1256 /* OK, if we didn't find a lower value we will have to locate one now - there must
1257 * be one in the interval [a=0,c].
1258 * The same thing is valid here, though: Don't spend dozens of iterations to find
1259 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1260 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1262 * I also have a safeguard for potentially really patological functions so we never
1263 * take more than 20 steps before we give up ...
1265 * If we already found a lower value we just skip this step and continue to the update.
1273 /* Select a new trial point.
1274 * If the derivatives at points a & c have different sign we interpolate to zero,
1275 * otherwise just do a bisection.
1277 if (gpa < 0 && gpc > 0)
1279 b = a + gpa*(a-c)/(gpc-gpa);
1286 /* safeguard if interpolation close to machine accuracy causes errors:
1287 * never go outside the interval
1289 if (b <= a || b >= c)
1294 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1296 /* Reload the old state */
1297 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1298 s_min, top, mdatoms, fr, vsite, constr,
1302 /* Take a trial step to this new point - new coords in s_b */
1303 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1304 constr, top, nrnb, wcycle, -1);
1307 /* Calculate energy for the trial step */
1308 evaluate_energy(fplog, cr,
1309 top_global, s_b, top,
1310 inputrec, nrnb, wcycle, gstat,
1311 vsite, constr, fcd, graph, mdatoms, fr,
1312 mu_tot, enerd, vir, pres, -1, FALSE);
1314 /* p does not change within a step, but since the domain decomposition
1315 * might change, we have to use cg_p of s_b here.
1320 for (i = 0; i < mdatoms->homenr; i++)
1322 for (m = 0; m < DIM; m++)
1324 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1327 /* Sum the gradient along the line across CPUs */
1330 gmx_sumd(1, &gpb, cr);
1335 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1336 s_a->epot, s_b->epot, s_c->epot, gpb);
1339 epot_repl = s_b->epot;
1341 /* Keep one of the intervals based on the value of the derivative at the new point */
1344 /* Replace c endpoint with b */
1345 swap_em_state(s_b, s_c);
1351 /* Replace a endpoint with b */
1352 swap_em_state(s_b, s_a);
1358 * Stop search as soon as we find a value smaller than the endpoints.
1359 * Never run more than 20 steps, no matter what.
1363 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1366 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1369 /* OK. We couldn't find a significantly lower energy.
1370 * If beta==0 this was steepest descent, and then we give up.
1371 * If not, set beta=0 and restart with steepest descent before quitting.
1381 /* Reset memory before giving up */
1387 /* Select min energy state of A & C, put the best in B.
1389 if (s_c->epot < s_a->epot)
1393 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1394 s_c->epot, s_a->epot);
1396 swap_em_state(s_b, s_c);
1404 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1405 s_a->epot, s_c->epot);
1407 swap_em_state(s_b, s_a);
1417 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1420 swap_em_state(s_b, s_c);
1425 /* new search direction */
1426 /* beta = 0 means forget all memory and restart with steepest descents. */
1427 if (nstcg && ((step % nstcg) == 0))
1433 /* s_min->fnorm cannot be zero, because then we would have converged
1437 /* Polak-Ribiere update.
1438 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1440 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1442 /* Limit beta to prevent oscillations */
1443 if (fabs(beta) > 5.0)
1449 /* update positions */
1450 swap_em_state(s_min, s_b);
1453 /* Print it if necessary */
1458 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1459 step, s_min->epot, s_min->fnorm/sqrt(state_global->natoms),
1460 s_min->fmax, s_min->a_fmax+1);
1462 /* Store the new (lower) energies */
1463 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1464 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1465 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1467 do_log = do_per_step(step, inputrec->nstlog);
1468 do_ene = do_per_step(step, inputrec->nstenergy);
1470 /* Prepare IMD energy record, if bIMD is TRUE. */
1471 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1475 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1477 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1478 do_log ? fplog : NULL, step, step, eprNORMAL,
1479 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1482 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1483 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
1485 IMD_send_positions(inputrec->imd);
1488 /* Stop when the maximum force lies below tolerance.
1489 * If we have reached machine precision, converged is already set to true.
1491 converged = converged || (s_min->fmax < inputrec->em_tol);
1493 } /* End of the loop */
1495 /* IMD cleanup, if bIMD is TRUE. */
1496 IMD_finalize(inputrec->bIMD, inputrec->imd);
1500 step--; /* we never took that last step in this case */
1503 if (s_min->fmax > inputrec->em_tol)
1507 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1508 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1515 /* If we printed energy and/or logfile last step (which was the last step)
1516 * we don't have to do it again, but otherwise print the final values.
1520 /* Write final value to log since we didn't do anything the last step */
1521 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1523 if (!do_ene || !do_log)
1525 /* Write final energy file entries */
1526 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1527 !do_log ? fplog : NULL, step, step, eprNORMAL,
1528 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1532 /* Print some stuff... */
1535 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1539 * For accurate normal mode calculation it is imperative that we
1540 * store the last conformation into the full precision binary trajectory.
1542 * However, we should only do it if we did NOT already write this step
1543 * above (which we did if do_x or do_f was true).
1545 do_x = !do_per_step(step, inputrec->nstxout);
1546 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1548 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1549 top_global, inputrec, step,
1550 s_min, state_global, f_global);
1552 fnormn = s_min->fnorm/sqrt(state_global->natoms);
1556 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1557 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1558 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1559 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1561 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1564 finish_em(cr, outf, walltime_accounting, wcycle);
1566 /* To print the actual number of steps we needed somewhere */
1567 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1570 } /* That's all folks */
1573 double do_lbfgs(FILE *fplog, t_commrec *cr,
1574 int nfile, const t_filenm fnm[],
1575 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
1576 int gmx_unused nstglobalcomm,
1577 gmx_vsite_t *vsite, gmx_constr_t constr,
1578 int gmx_unused stepout,
1579 t_inputrec *inputrec,
1580 gmx_mtop_t *top_global, t_fcdata *fcd,
1583 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1584 gmx_edsam_t gmx_unused ed,
1586 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1587 gmx_membed_t gmx_unused membed,
1588 real gmx_unused cpt_period, real gmx_unused max_hours,
1589 const char gmx_unused *deviceOptions,
1591 unsigned long gmx_unused Flags,
1592 gmx_walltime_accounting_t walltime_accounting)
1594 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1596 gmx_localtop_t *top;
1597 gmx_enerdata_t *enerd;
1599 gmx_global_stat_t gstat;
1602 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1603 double stepsize, gpa, gpb, gpc, tmp, minstep;
1604 real *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1605 real *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1606 real a, b, c, maxdelta, delta;
1607 real diag, Epot0, Epot, EpotA, EpotB, EpotC;
1608 real dgdx, dgdg, sq, yr, beta;
1610 gmx_bool converged, first;
1613 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1615 int start, end, number_steps;
1617 int i, k, m, n, nfmax, gf, step;
1624 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1629 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).");
1632 n = 3*state->natoms;
1633 nmaxcorr = inputrec->nbfgscorr;
1635 /* Allocate memory */
1636 /* Use pointers to real so we dont have to loop over both atoms and
1637 * dimensions all the time...
1638 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1639 * that point to the same memory.
1652 snew(rho, nmaxcorr);
1653 snew(alpha, nmaxcorr);
1656 for (i = 0; i < nmaxcorr; i++)
1662 for (i = 0; i < nmaxcorr; i++)
1671 init_em(fplog, LBFGS, cr, inputrec,
1672 state, top_global, &ems, &top, &f, &f_global,
1673 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1674 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
1675 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1676 * so we free some memory again.
1681 xx = (real *)state->x;
1685 end = mdatoms->homenr;
1687 /* Print to log file */
1688 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1690 do_log = do_ene = do_x = do_f = TRUE;
1692 /* Max number of steps */
1693 number_steps = inputrec->nsteps;
1695 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1697 for (i = start; i < end; i++)
1699 if (mdatoms->cFREEZE)
1701 gf = mdatoms->cFREEZE[i];
1703 for (m = 0; m < DIM; m++)
1705 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1710 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1714 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1719 construct_vsites(vsite, state->x, 1, NULL,
1720 top->idef.iparams, top->idef.il,
1721 fr->ePBC, fr->bMolPBC, cr, state->box);
1724 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1725 /* do_force always puts the charge groups in the box and shifts again
1726 * We do not unshift, so molecules are always whole
1731 evaluate_energy(fplog, cr,
1732 top_global, &ems, top,
1733 inputrec, nrnb, wcycle, gstat,
1734 vsite, constr, fcd, graph, mdatoms, fr,
1735 mu_tot, enerd, vir, pres, -1, TRUE);
1740 /* Copy stuff to the energy bin for easy printing etc. */
1741 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1742 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
1743 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1745 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
1746 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1747 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1751 /* This is the starting energy */
1752 Epot = enerd->term[F_EPOT];
1758 /* Set the initial step.
1759 * since it will be multiplied by the non-normalized search direction
1760 * vector (force vector the first time), we scale it by the
1761 * norm of the force.
1766 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1767 fprintf(stderr, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1768 fprintf(stderr, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1769 fprintf(stderr, "\n");
1770 /* and copy to the log file too... */
1771 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1772 fprintf(fplog, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1773 fprintf(fplog, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1774 fprintf(fplog, "\n");
1778 for (i = 0; i < n; i++)
1782 dx[point][i] = ff[i]; /* Initial search direction */
1790 stepsize = 1.0/fnorm;
1792 /* Start the loop over BFGS steps.
1793 * Each successful step is counted, and we continue until
1794 * we either converge or reach the max number of steps.
1799 /* Set the gradient from the force */
1801 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1804 /* Write coordinates if necessary */
1805 do_x = do_per_step(step, inputrec->nstxout);
1806 do_f = do_per_step(step, inputrec->nstfout);
1811 mdof_flags |= MDOF_X;
1816 mdof_flags |= MDOF_F;
1821 mdof_flags |= MDOF_IMD;
1824 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1825 top_global, step, (real)step, state, state, f, f);
1827 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1829 /* pointer to current direction - point=0 first time here */
1832 /* calculate line gradient */
1833 for (gpa = 0, i = 0; i < n; i++)
1838 /* Calculate minimum allowed stepsize, before the average (norm)
1839 * relative change in coordinate is smaller than precision
1841 for (minstep = 0, i = 0; i < n; i++)
1851 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1853 if (stepsize < minstep)
1859 /* Store old forces and coordinates */
1860 for (i = 0; i < n; i++)
1869 for (i = 0; i < n; i++)
1874 /* Take a step downhill.
1875 * In theory, we should minimize the function along this direction.
1876 * That is quite possible, but it turns out to take 5-10 function evaluations
1877 * for each line. However, we dont really need to find the exact minimum -
1878 * it is much better to start a new BFGS step in a modified direction as soon
1879 * as we are close to it. This will save a lot of energy evaluations.
1881 * In practice, we just try to take a single step.
1882 * If it worked (i.e. lowered the energy), we increase the stepsize but
1883 * the continue straight to the next BFGS step without trying to find any minimum.
1884 * If it didn't work (higher energy), there must be a minimum somewhere between
1885 * the old position and the new one.
1887 * Due to the finite numerical accuracy, it turns out that it is a good idea
1888 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1889 * This leads to lower final energies in the tests I've done. / Erik
1894 c = a + stepsize; /* reference position along line is zero */
1896 /* Check stepsize first. We do not allow displacements
1897 * larger than emstep.
1903 for (i = 0; i < n; i++)
1906 if (delta > maxdelta)
1911 if (maxdelta > inputrec->em_stepsize)
1916 while (maxdelta > inputrec->em_stepsize);
1918 /* Take a trial step */
1919 for (i = 0; i < n; i++)
1921 xc[i] = lastx[i] + c*s[i];
1925 /* Calculate energy for the trial step */
1926 ems.s.x = (rvec *)xc;
1928 evaluate_energy(fplog, cr,
1929 top_global, &ems, top,
1930 inputrec, nrnb, wcycle, gstat,
1931 vsite, constr, fcd, graph, mdatoms, fr,
1932 mu_tot, enerd, vir, pres, step, FALSE);
1935 /* Calc derivative along line */
1936 for (gpc = 0, i = 0; i < n; i++)
1938 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1940 /* Sum the gradient along the line across CPUs */
1943 gmx_sumd(1, &gpc, cr);
1946 /* This is the max amount of increase in energy we tolerate */
1947 tmp = sqrt(GMX_REAL_EPS)*fabs(EpotA);
1949 /* Accept the step if the energy is lower, or if it is not significantly higher
1950 * and the line derivative is still negative.
1952 if (EpotC < EpotA || (gpc < 0 && EpotC < (EpotA+tmp)))
1955 /* Great, we found a better energy. Increase step for next iteration
1956 * if we are still going down, decrease it otherwise
1960 stepsize *= 1.618034; /* The golden section */
1964 stepsize *= 0.618034; /* 1/golden section */
1969 /* New energy is the same or higher. We will have to do some work
1970 * to find a smaller value in the interval. Take smaller step next time!
1973 stepsize *= 0.618034;
1976 /* OK, if we didn't find a lower value we will have to locate one now - there must
1977 * be one in the interval [a=0,c].
1978 * The same thing is valid here, though: Don't spend dozens of iterations to find
1979 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1980 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1982 * I also have a safeguard for potentially really patological functions so we never
1983 * take more than 20 steps before we give up ...
1985 * If we already found a lower value we just skip this step and continue to the update.
1994 /* Select a new trial point.
1995 * If the derivatives at points a & c have different sign we interpolate to zero,
1996 * otherwise just do a bisection.
1999 if (gpa < 0 && gpc > 0)
2001 b = a + gpa*(a-c)/(gpc-gpa);
2008 /* safeguard if interpolation close to machine accuracy causes errors:
2009 * never go outside the interval
2011 if (b <= a || b >= c)
2016 /* Take a trial step */
2017 for (i = 0; i < n; i++)
2019 xb[i] = lastx[i] + b*s[i];
2023 /* Calculate energy for the trial step */
2024 ems.s.x = (rvec *)xb;
2026 evaluate_energy(fplog, cr,
2027 top_global, &ems, top,
2028 inputrec, nrnb, wcycle, gstat,
2029 vsite, constr, fcd, graph, mdatoms, fr,
2030 mu_tot, enerd, vir, pres, step, FALSE);
2035 for (gpb = 0, i = 0; i < n; i++)
2037 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2040 /* Sum the gradient along the line across CPUs */
2043 gmx_sumd(1, &gpb, cr);
2046 /* Keep one of the intervals based on the value of the derivative at the new point */
2049 /* Replace c endpoint with b */
2053 /* swap coord pointers b/c */
2063 /* Replace a endpoint with b */
2067 /* swap coord pointers a/b */
2077 * Stop search as soon as we find a value smaller than the endpoints,
2078 * or if the tolerance is below machine precision.
2079 * Never run more than 20 steps, no matter what.
2083 while ((EpotB > EpotA || EpotB > EpotC) && (nminstep < 20));
2085 if (fabs(EpotB-Epot0) < GMX_REAL_EPS || nminstep >= 20)
2087 /* OK. We couldn't find a significantly lower energy.
2088 * If ncorr==0 this was steepest descent, and then we give up.
2089 * If not, reset memory to restart as steepest descent before quitting.
2101 /* Search in gradient direction */
2102 for (i = 0; i < n; i++)
2104 dx[point][i] = ff[i];
2106 /* Reset stepsize */
2107 stepsize = 1.0/fnorm;
2112 /* Select min energy state of A & C, put the best in xx/ff/Epot
2118 for (i = 0; i < n; i++)
2129 for (i = 0; i < n; i++)
2143 for (i = 0; i < n; i++)
2151 /* Update the memory information, and calculate a new
2152 * approximation of the inverse hessian
2155 /* Have new data in Epot, xx, ff */
2156 if (ncorr < nmaxcorr)
2161 for (i = 0; i < n; i++)
2163 dg[point][i] = lastf[i]-ff[i];
2164 dx[point][i] *= stepsize;
2169 for (i = 0; i < n; i++)
2171 dgdg += dg[point][i]*dg[point][i];
2172 dgdx += dg[point][i]*dx[point][i];
2177 rho[point] = 1.0/dgdx;
2180 if (point >= nmaxcorr)
2186 for (i = 0; i < n; i++)
2193 /* Recursive update. First go back over the memory points */
2194 for (k = 0; k < ncorr; k++)
2203 for (i = 0; i < n; i++)
2205 sq += dx[cp][i]*p[i];
2208 alpha[cp] = rho[cp]*sq;
2210 for (i = 0; i < n; i++)
2212 p[i] -= alpha[cp]*dg[cp][i];
2216 for (i = 0; i < n; i++)
2221 /* And then go forward again */
2222 for (k = 0; k < ncorr; k++)
2225 for (i = 0; i < n; i++)
2227 yr += p[i]*dg[cp][i];
2231 beta = alpha[cp]-beta;
2233 for (i = 0; i < n; i++)
2235 p[i] += beta*dx[cp][i];
2245 for (i = 0; i < n; i++)
2249 dx[point][i] = p[i];
2259 /* Test whether the convergence criterion is met */
2260 get_f_norm_max(cr, &(inputrec->opts), mdatoms, f, &fnorm, &fmax, &nfmax);
2262 /* Print it if necessary */
2267 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2268 step, Epot, fnorm/sqrt(state->natoms), fmax, nfmax+1);
2270 /* Store the new (lower) energies */
2271 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2272 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
2273 NULL, NULL, vir, pres, NULL, mu_tot, constr);
2274 do_log = do_per_step(step, inputrec->nstlog);
2275 do_ene = do_per_step(step, inputrec->nstenergy);
2278 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2280 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2281 do_log ? fplog : NULL, step, step, eprNORMAL,
2282 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2285 /* Send x and E to IMD client, if bIMD is TRUE. */
2286 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state->box, state->x, inputrec, 0, wcycle) && MASTER(cr))
2288 IMD_send_positions(inputrec->imd);
2291 /* Stop when the maximum force lies below tolerance.
2292 * If we have reached machine precision, converged is already set to true.
2295 converged = converged || (fmax < inputrec->em_tol);
2297 } /* End of the loop */
2299 /* IMD cleanup, if bIMD is TRUE. */
2300 IMD_finalize(inputrec->bIMD, inputrec->imd);
2304 step--; /* we never took that last step in this case */
2307 if (fmax > inputrec->em_tol)
2311 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2312 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2317 /* If we printed energy and/or logfile last step (which was the last step)
2318 * we don't have to do it again, but otherwise print the final values.
2320 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2322 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2324 if (!do_ene || !do_log) /* Write final energy file entries */
2326 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2327 !do_log ? fplog : NULL, step, step, eprNORMAL,
2328 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2331 /* Print some stuff... */
2334 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2338 * For accurate normal mode calculation it is imperative that we
2339 * store the last conformation into the full precision binary trajectory.
2341 * However, we should only do it if we did NOT already write this step
2342 * above (which we did if do_x or do_f was true).
2344 do_x = !do_per_step(step, inputrec->nstxout);
2345 do_f = !do_per_step(step, inputrec->nstfout);
2346 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2347 top_global, inputrec, step,
2352 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2353 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2354 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2355 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2357 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2360 finish_em(cr, outf, walltime_accounting, wcycle);
2362 /* To print the actual number of steps we needed somewhere */
2363 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2366 } /* That's all folks */
2369 double do_steep(FILE *fplog, t_commrec *cr,
2370 int nfile, const t_filenm fnm[],
2371 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2372 int gmx_unused nstglobalcomm,
2373 gmx_vsite_t *vsite, gmx_constr_t constr,
2374 int gmx_unused stepout,
2375 t_inputrec *inputrec,
2376 gmx_mtop_t *top_global, t_fcdata *fcd,
2377 t_state *state_global,
2379 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2380 gmx_edsam_t gmx_unused ed,
2382 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2383 gmx_membed_t gmx_unused membed,
2384 real gmx_unused cpt_period, real gmx_unused max_hours,
2385 const char gmx_unused *deviceOptions,
2387 unsigned long gmx_unused Flags,
2388 gmx_walltime_accounting_t walltime_accounting)
2390 const char *SD = "Steepest Descents";
2391 em_state_t *s_min, *s_try;
2393 gmx_localtop_t *top;
2394 gmx_enerdata_t *enerd;
2396 gmx_global_stat_t gstat;
2398 real stepsize, constepsize;
2402 gmx_bool bDone, bAbort, do_x, do_f;
2407 int steps_accepted = 0;
2411 s_min = init_em_state();
2412 s_try = init_em_state();
2414 /* Init em and store the local state in s_try */
2415 init_em(fplog, SD, cr, inputrec,
2416 state_global, top_global, s_try, &top, &f, &f_global,
2417 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2418 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
2420 /* Print to log file */
2421 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2423 /* Set variables for stepsize (in nm). This is the largest
2424 * step that we are going to make in any direction.
2426 ustep = inputrec->em_stepsize;
2429 /* Max number of steps */
2430 nsteps = inputrec->nsteps;
2434 /* Print to the screen */
2435 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2439 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2442 /**** HERE STARTS THE LOOP ****
2443 * count is the counter for the number of steps
2444 * bDone will be TRUE when the minimization has converged
2445 * bAbort will be TRUE when nsteps steps have been performed or when
2446 * the stepsize becomes smaller than is reasonable for machine precision
2451 while (!bDone && !bAbort)
2453 bAbort = (nsteps >= 0) && (count == nsteps);
2455 /* set new coordinates, except for first step */
2458 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2459 s_min, stepsize, s_min->f, s_try,
2460 constr, top, nrnb, wcycle, count);
2463 evaluate_energy(fplog, cr,
2464 top_global, s_try, top,
2465 inputrec, nrnb, wcycle, gstat,
2466 vsite, constr, fcd, graph, mdatoms, fr,
2467 mu_tot, enerd, vir, pres, count, count == 0);
2471 print_ebin_header(fplog, count, count, s_try->s.lambda[efptFEP]);
2476 s_min->epot = s_try->epot + 1;
2479 /* Print it if necessary */
2484 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2485 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2486 (s_try->epot < s_min->epot) ? '\n' : '\r');
2489 if (s_try->epot < s_min->epot)
2491 /* Store the new (lower) energies */
2492 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2493 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2494 s_try->s.box, NULL, NULL, vir, pres, NULL, mu_tot, constr);
2496 /* Prepare IMD energy record, if bIMD is TRUE. */
2497 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2499 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2500 do_per_step(steps_accepted, inputrec->nstdisreout),
2501 do_per_step(steps_accepted, inputrec->nstorireout),
2502 fplog, count, count, eprNORMAL, TRUE,
2503 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2508 /* Now if the new energy is smaller than the previous...
2509 * or if this is the first step!
2510 * or if we did random steps!
2513 if ( (count == 0) || (s_try->epot < s_min->epot) )
2517 /* Test whether the convergence criterion is met... */
2518 bDone = (s_try->fmax < inputrec->em_tol);
2520 /* Copy the arrays for force, positions and energy */
2521 /* The 'Min' array always holds the coords and forces of the minimal
2523 swap_em_state(s_min, s_try);
2529 /* Write to trn, if necessary */
2530 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2531 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2532 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
2533 top_global, inputrec, count,
2534 s_min, state_global, f_global);
2538 /* If energy is not smaller make the step smaller... */
2541 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2543 /* Reload the old state */
2544 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2545 s_min, top, mdatoms, fr, vsite, constr,
2550 /* Determine new step */
2551 stepsize = ustep/s_min->fmax;
2553 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2555 if (count == nsteps || ustep < 1e-12)
2557 if (count == nsteps || ustep < 1e-6)
2562 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2563 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2568 /* Send IMD energies and positions, if bIMD is TRUE. */
2569 if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
2571 IMD_send_positions(inputrec->imd);
2575 } /* End of the loop */
2577 /* IMD cleanup, if bIMD is TRUE. */
2578 IMD_finalize(inputrec->bIMD, inputrec->imd);
2580 /* Print some data... */
2583 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2585 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2586 top_global, inputrec, count,
2587 s_min, state_global, f_global);
2589 fnormn = s_min->fnorm/sqrt(state_global->natoms);
2593 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2594 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2595 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2596 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2599 finish_em(cr, outf, walltime_accounting, wcycle);
2601 /* To print the actual number of steps we needed somewhere */
2602 inputrec->nsteps = count;
2604 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2607 } /* That's all folks */
2610 double do_nm(FILE *fplog, t_commrec *cr,
2611 int nfile, const t_filenm fnm[],
2612 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2613 int gmx_unused nstglobalcomm,
2614 gmx_vsite_t *vsite, gmx_constr_t constr,
2615 int gmx_unused stepout,
2616 t_inputrec *inputrec,
2617 gmx_mtop_t *top_global, t_fcdata *fcd,
2618 t_state *state_global,
2620 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2621 gmx_edsam_t gmx_unused ed,
2623 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2624 gmx_membed_t gmx_unused membed,
2625 real gmx_unused cpt_period, real gmx_unused max_hours,
2626 const char gmx_unused *deviceOptions,
2628 unsigned long gmx_unused Flags,
2629 gmx_walltime_accounting_t walltime_accounting)
2631 const char *NM = "Normal Mode Analysis";
2633 int natoms, atom, d;
2636 gmx_localtop_t *top;
2637 gmx_enerdata_t *enerd;
2639 gmx_global_stat_t gstat;
2641 real t, t0, lambda, lam0;
2646 gmx_bool bSparse; /* use sparse matrix storage format */
2648 gmx_sparsematrix_t * sparse_matrix = NULL;
2649 real * full_matrix = NULL;
2650 em_state_t * state_work;
2652 /* added with respect to mdrun */
2653 int i, j, k, row, col;
2654 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2660 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2663 state_work = init_em_state();
2665 /* Init em and store the local state in state_minimum */
2666 init_em(fplog, NM, cr, inputrec,
2667 state_global, top_global, state_work, &top,
2669 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2670 nfile, fnm, &outf, NULL, imdport, Flags, wcycle);
2672 natoms = top_global->natoms;
2680 "NOTE: This version of Gromacs has been compiled in single precision,\n"
2681 " which MIGHT not be accurate enough for normal mode analysis.\n"
2682 " Gromacs now uses sparse matrix storage, so the memory requirements\n"
2683 " are fairly modest even if you recompile in double precision.\n\n");
2687 /* Check if we can/should use sparse storage format.
2689 * Sparse format is only useful when the Hessian itself is sparse, which it
2690 * will be when we use a cutoff.
2691 * For small systems (n<1000) it is easier to always use full matrix format, though.
2693 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2695 md_print_info(cr, fplog, "Non-cutoff electrostatics used, forcing full Hessian format.\n");
2698 else if (top_global->natoms < 1000)
2700 md_print_info(cr, fplog, "Small system size (N=%d), using full Hessian format.\n", top_global->natoms);
2705 md_print_info(cr, fplog, "Using compressed symmetric sparse Hessian format.\n");
2711 sz = DIM*top_global->natoms;
2713 fprintf(stderr, "Allocating Hessian memory...\n\n");
2717 sparse_matrix = gmx_sparsematrix_init(sz);
2718 sparse_matrix->compressed_symmetric = TRUE;
2722 snew(full_matrix, sz*sz);
2726 /* Initial values */
2727 t0 = inputrec->init_t;
2728 lam0 = inputrec->fepvals->init_lambda;
2736 /* Write start time and temperature */
2737 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2739 /* fudge nr of steps to nr of atoms */
2740 inputrec->nsteps = natoms*2;
2744 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2745 *(top_global->name), (int)inputrec->nsteps);
2748 nnodes = cr->nnodes;
2750 /* Make evaluate_energy do a single node force calculation */
2752 evaluate_energy(fplog, cr,
2753 top_global, state_work, top,
2754 inputrec, nrnb, wcycle, gstat,
2755 vsite, constr, fcd, graph, mdatoms, fr,
2756 mu_tot, enerd, vir, pres, -1, TRUE);
2757 cr->nnodes = nnodes;
2759 /* if forces are not small, warn user */
2760 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, state_work);
2762 md_print_info(cr, fplog, "Maximum force:%12.5e\n", state_work->fmax);
2763 if (state_work->fmax > 1.0e-3)
2765 md_print_info(cr, fplog,
2766 "The force is probably not small enough to "
2767 "ensure that you are at a minimum.\n"
2768 "Be aware that negative eigenvalues may occur\n"
2769 "when the resulting matrix is diagonalized.\n\n");
2772 /***********************************************************
2774 * Loop over all pairs in matrix
2776 * do_force called twice. Once with positive and
2777 * once with negative displacement
2779 ************************************************************/
2781 /* Steps are divided one by one over the nodes */
2782 for (atom = cr->nodeid; atom < natoms; atom += nnodes)
2785 for (d = 0; d < DIM; d++)
2787 x_min = state_work->s.x[atom][d];
2789 state_work->s.x[atom][d] = x_min - der_range;
2791 /* Make evaluate_energy do a single node force calculation */
2793 evaluate_energy(fplog, cr,
2794 top_global, state_work, top,
2795 inputrec, nrnb, wcycle, gstat,
2796 vsite, constr, fcd, graph, mdatoms, fr,
2797 mu_tot, enerd, vir, pres, atom*2, FALSE);
2799 for (i = 0; i < natoms; i++)
2801 copy_rvec(state_work->f[i], fneg[i]);
2804 state_work->s.x[atom][d] = x_min + der_range;
2806 evaluate_energy(fplog, cr,
2807 top_global, state_work, top,
2808 inputrec, nrnb, wcycle, gstat,
2809 vsite, constr, fcd, graph, mdatoms, fr,
2810 mu_tot, enerd, vir, pres, atom*2+1, FALSE);
2811 cr->nnodes = nnodes;
2813 /* x is restored to original */
2814 state_work->s.x[atom][d] = x_min;
2816 for (j = 0; j < natoms; j++)
2818 for (k = 0; (k < DIM); k++)
2821 -(state_work->f[j][k] - fneg[j][k])/(2*der_range);
2829 #define mpi_type MPI_DOUBLE
2831 #define mpi_type MPI_FLOAT
2833 MPI_Send(dfdx[0], natoms*DIM, mpi_type, MASTERNODE(cr), cr->nodeid,
2834 cr->mpi_comm_mygroup);
2839 for (node = 0; (node < nnodes && atom+node < natoms); node++)
2845 MPI_Recv(dfdx[0], natoms*DIM, mpi_type, node, node,
2846 cr->mpi_comm_mygroup, &stat);
2851 row = (atom + node)*DIM + d;
2853 for (j = 0; j < natoms; j++)
2855 for (k = 0; k < DIM; k++)
2861 if (col >= row && dfdx[j][k] != 0.0)
2863 gmx_sparsematrix_increment_value(sparse_matrix,
2864 row, col, dfdx[j][k]);
2869 full_matrix[row*sz+col] = dfdx[j][k];
2876 if (bVerbose && fplog)
2881 /* write progress */
2882 if (MASTER(cr) && bVerbose)
2884 fprintf(stderr, "\rFinished step %d out of %d",
2885 min(atom+nnodes, natoms), natoms);
2892 fprintf(stderr, "\n\nWriting Hessian...\n");
2893 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2896 finish_em(cr, outf, walltime_accounting, wcycle);
2898 walltime_accounting_set_nsteps_done(walltime_accounting, natoms*2);