2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
59 #include "md_support.h"
64 #include "gromacs/topology/mtop_util.h"
67 #include "gmx_omp_nthreads.h"
68 #include "md_logging.h"
70 #include "gromacs/fileio/confio.h"
71 #include "gromacs/fileio/mtxio.h"
72 #include "gromacs/fileio/trajectory_writing.h"
73 #include "gromacs/imd/imd.h"
74 #include "gromacs/legacyheaders/types/commrec.h"
75 #include "gromacs/linearalgebra/sparsematrix.h"
76 #include "gromacs/math/vec.h"
77 #include "gromacs/pbcutil/mshift.h"
78 #include "gromacs/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, 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, 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(fplog, inputrec, fr, count, 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, 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 if (fr->bSepDVDL && fplog)
796 gmx_print_sepdvdl(fplog, "Constraints", t, dvdl_constr);
798 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
799 m_add(force_vir, shake_vir, vir);
800 wallcycle_stop(wcycle, ewcCONSTR);
804 copy_mat(force_vir, vir);
808 enerd->term[F_PRES] =
809 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
811 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
813 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
815 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
819 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
821 em_state_t *s_min, em_state_t *s_b)
825 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
827 unsigned char *grpnrFREEZE;
831 fprintf(debug, "Doing reorder_partsum\n");
837 cgs_gl = dd_charge_groups_global(cr->dd);
838 index = cgs_gl->index;
840 /* Collect fm in a global vector fmg.
841 * This conflicts with the spirit of domain decomposition,
842 * but to fully optimize this a much more complicated algorithm is required.
844 snew(fmg, mtop->natoms);
846 ncg = s_min->s.ncg_gl;
847 cg_gl = s_min->s.cg_gl;
849 for (c = 0; c < ncg; c++)
854 for (a = a0; a < a1; a++)
856 copy_rvec(fm[i], fmg[a]);
860 gmx_sum(mtop->natoms*3, fmg[0], cr);
862 /* Now we will determine the part of the sum for the cgs in state s_b */
864 cg_gl = s_b->s.cg_gl;
868 grpnrFREEZE = mtop->groups.grpnr[egcFREEZE];
869 for (c = 0; c < ncg; c++)
874 for (a = a0; a < a1; a++)
876 if (mdatoms->cFREEZE && grpnrFREEZE)
880 for (m = 0; m < DIM; m++)
882 if (!opts->nFreeze[gf][m])
884 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
896 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
898 em_state_t *s_min, em_state_t *s_b)
904 /* This is just the classical Polak-Ribiere calculation of beta;
905 * it looks a bit complicated since we take freeze groups into account,
906 * and might have to sum it in parallel runs.
909 if (!DOMAINDECOMP(cr) ||
910 (s_min->s.ddp_count == cr->dd->ddp_count &&
911 s_b->s.ddp_count == cr->dd->ddp_count))
917 /* This part of code can be incorrect with DD,
918 * since the atom ordering in s_b and s_min might differ.
920 for (i = 0; i < mdatoms->homenr; i++)
922 if (mdatoms->cFREEZE)
924 gf = mdatoms->cFREEZE[i];
926 for (m = 0; m < DIM; m++)
928 if (!opts->nFreeze[gf][m])
930 sum += (fb[i][m] - fm[i][m])*fb[i][m];
937 /* We need to reorder cgs while summing */
938 sum = reorder_partsum(cr, opts, mdatoms, mtop, s_min, s_b);
942 gmx_sumd(1, &sum, cr);
945 return sum/sqr(s_min->fnorm);
948 double do_cg(FILE *fplog, t_commrec *cr,
949 int nfile, const t_filenm fnm[],
950 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
951 int gmx_unused nstglobalcomm,
952 gmx_vsite_t *vsite, gmx_constr_t constr,
953 int gmx_unused stepout,
954 t_inputrec *inputrec,
955 gmx_mtop_t *top_global, t_fcdata *fcd,
956 t_state *state_global,
958 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
959 gmx_edsam_t gmx_unused ed,
961 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
962 gmx_membed_t gmx_unused membed,
963 real gmx_unused cpt_period, real gmx_unused max_hours,
964 const char gmx_unused *deviceOptions,
966 unsigned long gmx_unused Flags,
967 gmx_walltime_accounting_t walltime_accounting)
969 const char *CG = "Polak-Ribiere Conjugate Gradients";
971 em_state_t *s_min, *s_a, *s_b, *s_c;
973 gmx_enerdata_t *enerd;
975 gmx_global_stat_t gstat;
977 rvec *f_global, *p, *sf, *sfm;
978 double gpa, gpb, gpc, tmp, sum[2], minstep;
981 real a, b, c, beta = 0.0;
985 gmx_bool converged, foundlower;
987 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
989 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
991 int i, m, gf, step, nminstep;
996 s_min = init_em_state();
997 s_a = init_em_state();
998 s_b = init_em_state();
999 s_c = init_em_state();
1001 /* Init em and store the local state in s_min */
1002 init_em(fplog, CG, cr, inputrec,
1003 state_global, top_global, s_min, &top, &f, &f_global,
1004 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1005 nfile, fnm, &outf, &mdebin, imdport, Flags);
1007 /* Print to log file */
1008 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1010 /* Max number of steps */
1011 number_steps = inputrec->nsteps;
1015 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1019 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1022 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1023 /* do_force always puts the charge groups in the box and shifts again
1024 * We do not unshift, so molecules are always whole in congrad.c
1026 evaluate_energy(fplog, cr,
1027 top_global, s_min, top,
1028 inputrec, nrnb, wcycle, gstat,
1029 vsite, constr, fcd, graph, mdatoms, fr,
1030 mu_tot, enerd, vir, pres, -1, TRUE);
1035 /* Copy stuff to the energy bin for easy printing etc. */
1036 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1037 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1038 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1040 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1041 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1042 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1046 /* Estimate/guess the initial stepsize */
1047 stepsize = inputrec->em_stepsize/s_min->fnorm;
1051 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1052 s_min->fmax, s_min->a_fmax+1);
1053 fprintf(stderr, " F-Norm = %12.5e\n",
1054 s_min->fnorm/sqrt(state_global->natoms));
1055 fprintf(stderr, "\n");
1056 /* and copy to the log file too... */
1057 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1058 s_min->fmax, s_min->a_fmax+1);
1059 fprintf(fplog, " F-Norm = %12.5e\n",
1060 s_min->fnorm/sqrt(state_global->natoms));
1061 fprintf(fplog, "\n");
1063 /* Start the loop over CG steps.
1064 * Each successful step is counted, and we continue until
1065 * we either converge or reach the max number of steps.
1068 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1071 /* start taking steps in a new direction
1072 * First time we enter the routine, beta=0, and the direction is
1073 * simply the negative gradient.
1076 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1081 for (i = 0; i < mdatoms->homenr; i++)
1083 if (mdatoms->cFREEZE)
1085 gf = mdatoms->cFREEZE[i];
1087 for (m = 0; m < DIM; m++)
1089 if (!inputrec->opts.nFreeze[gf][m])
1091 p[i][m] = sf[i][m] + beta*p[i][m];
1092 gpa -= p[i][m]*sf[i][m];
1093 /* f is negative gradient, thus the sign */
1102 /* Sum the gradient along the line across CPUs */
1105 gmx_sumd(1, &gpa, cr);
1108 /* Calculate the norm of the search vector */
1109 get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1111 /* Just in case stepsize reaches zero due to numerical precision... */
1114 stepsize = inputrec->em_stepsize/pnorm;
1118 * Double check the value of the derivative in the search direction.
1119 * If it is positive it must be due to the old information in the
1120 * CG formula, so just remove that and start over with beta=0.
1121 * This corresponds to a steepest descent step.
1126 step--; /* Don't count this step since we are restarting */
1127 continue; /* Go back to the beginning of the big for-loop */
1130 /* Calculate minimum allowed stepsize, before the average (norm)
1131 * relative change in coordinate is smaller than precision
1134 for (i = 0; i < mdatoms->homenr; i++)
1136 for (m = 0; m < DIM; m++)
1138 tmp = fabs(s_min->s.x[i][m]);
1147 /* Add up from all CPUs */
1150 gmx_sumd(1, &minstep, cr);
1153 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1155 if (stepsize < minstep)
1161 /* Write coordinates if necessary */
1162 do_x = do_per_step(step, inputrec->nstxout);
1163 do_f = do_per_step(step, inputrec->nstfout);
1165 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1166 top_global, inputrec, step,
1167 s_min, state_global, f_global);
1169 /* Take a step downhill.
1170 * In theory, we should minimize the function along this direction.
1171 * That is quite possible, but it turns out to take 5-10 function evaluations
1172 * for each line. However, we dont really need to find the exact minimum -
1173 * it is much better to start a new CG step in a modified direction as soon
1174 * as we are close to it. This will save a lot of energy evaluations.
1176 * In practice, we just try to take a single step.
1177 * If it worked (i.e. lowered the energy), we increase the stepsize but
1178 * the continue straight to the next CG step without trying to find any minimum.
1179 * If it didn't work (higher energy), there must be a minimum somewhere between
1180 * the old position and the new one.
1182 * Due to the finite numerical accuracy, it turns out that it is a good idea
1183 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1184 * This leads to lower final energies in the tests I've done. / Erik
1186 s_a->epot = s_min->epot;
1188 c = a + stepsize; /* reference position along line is zero */
1190 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1192 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1193 s_min, top, mdatoms, fr, vsite, constr,
1197 /* Take a trial step (new coords in s_c) */
1198 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1199 constr, top, nrnb, wcycle, -1);
1202 /* Calculate energy for the trial step */
1203 evaluate_energy(fplog, cr,
1204 top_global, s_c, top,
1205 inputrec, nrnb, wcycle, gstat,
1206 vsite, constr, fcd, graph, mdatoms, fr,
1207 mu_tot, enerd, vir, pres, -1, FALSE);
1209 /* Calc derivative along line */
1213 for (i = 0; i < mdatoms->homenr; i++)
1215 for (m = 0; m < DIM; m++)
1217 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1220 /* Sum the gradient along the line across CPUs */
1223 gmx_sumd(1, &gpc, cr);
1226 /* This is the max amount of increase in energy we tolerate */
1227 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1229 /* Accept the step if the energy is lower, or if it is not significantly higher
1230 * and the line derivative is still negative.
1232 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1235 /* Great, we found a better energy. Increase step for next iteration
1236 * if we are still going down, decrease it otherwise
1240 stepsize *= 1.618034; /* The golden section */
1244 stepsize *= 0.618034; /* 1/golden section */
1249 /* New energy is the same or higher. We will have to do some work
1250 * to find a smaller value in the interval. Take smaller step next time!
1253 stepsize *= 0.618034;
1259 /* OK, if we didn't find a lower value we will have to locate one now - there must
1260 * be one in the interval [a=0,c].
1261 * The same thing is valid here, though: Don't spend dozens of iterations to find
1262 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1263 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1265 * I also have a safeguard for potentially really patological functions so we never
1266 * take more than 20 steps before we give up ...
1268 * If we already found a lower value we just skip this step and continue to the update.
1276 /* Select a new trial point.
1277 * If the derivatives at points a & c have different sign we interpolate to zero,
1278 * otherwise just do a bisection.
1280 if (gpa < 0 && gpc > 0)
1282 b = a + gpa*(a-c)/(gpc-gpa);
1289 /* safeguard if interpolation close to machine accuracy causes errors:
1290 * never go outside the interval
1292 if (b <= a || b >= c)
1297 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1299 /* Reload the old state */
1300 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1301 s_min, top, mdatoms, fr, vsite, constr,
1305 /* Take a trial step to this new point - new coords in s_b */
1306 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1307 constr, top, nrnb, wcycle, -1);
1310 /* Calculate energy for the trial step */
1311 evaluate_energy(fplog, cr,
1312 top_global, s_b, top,
1313 inputrec, nrnb, wcycle, gstat,
1314 vsite, constr, fcd, graph, mdatoms, fr,
1315 mu_tot, enerd, vir, pres, -1, FALSE);
1317 /* p does not change within a step, but since the domain decomposition
1318 * might change, we have to use cg_p of s_b here.
1323 for (i = 0; i < mdatoms->homenr; i++)
1325 for (m = 0; m < DIM; m++)
1327 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1330 /* Sum the gradient along the line across CPUs */
1333 gmx_sumd(1, &gpb, cr);
1338 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1339 s_a->epot, s_b->epot, s_c->epot, gpb);
1342 epot_repl = s_b->epot;
1344 /* Keep one of the intervals based on the value of the derivative at the new point */
1347 /* Replace c endpoint with b */
1348 swap_em_state(s_b, s_c);
1354 /* Replace a endpoint with b */
1355 swap_em_state(s_b, s_a);
1361 * Stop search as soon as we find a value smaller than the endpoints.
1362 * Never run more than 20 steps, no matter what.
1366 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1369 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1372 /* OK. We couldn't find a significantly lower energy.
1373 * If beta==0 this was steepest descent, and then we give up.
1374 * If not, set beta=0 and restart with steepest descent before quitting.
1384 /* Reset memory before giving up */
1390 /* Select min energy state of A & C, put the best in B.
1392 if (s_c->epot < s_a->epot)
1396 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1397 s_c->epot, s_a->epot);
1399 swap_em_state(s_b, s_c);
1407 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1408 s_a->epot, s_c->epot);
1410 swap_em_state(s_b, s_a);
1420 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1423 swap_em_state(s_b, s_c);
1428 /* new search direction */
1429 /* beta = 0 means forget all memory and restart with steepest descents. */
1430 if (nstcg && ((step % nstcg) == 0))
1436 /* s_min->fnorm cannot be zero, because then we would have converged
1440 /* Polak-Ribiere update.
1441 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1443 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1445 /* Limit beta to prevent oscillations */
1446 if (fabs(beta) > 5.0)
1452 /* update positions */
1453 swap_em_state(s_min, s_b);
1456 /* Print it if necessary */
1461 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1462 step, s_min->epot, s_min->fnorm/sqrt(state_global->natoms),
1463 s_min->fmax, s_min->a_fmax+1);
1465 /* Store the new (lower) energies */
1466 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1467 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1468 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1470 do_log = do_per_step(step, inputrec->nstlog);
1471 do_ene = do_per_step(step, inputrec->nstenergy);
1473 /* Prepare IMD energy record, if bIMD is TRUE. */
1474 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1478 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1480 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1481 do_log ? fplog : NULL, step, step, eprNORMAL,
1482 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1485 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1486 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
1488 IMD_send_positions(inputrec->imd);
1491 /* Stop when the maximum force lies below tolerance.
1492 * If we have reached machine precision, converged is already set to true.
1494 converged = converged || (s_min->fmax < inputrec->em_tol);
1496 } /* End of the loop */
1498 /* IMD cleanup, if bIMD is TRUE. */
1499 IMD_finalize(inputrec->bIMD, inputrec->imd);
1503 step--; /* we never took that last step in this case */
1506 if (s_min->fmax > inputrec->em_tol)
1510 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1511 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1518 /* If we printed energy and/or logfile last step (which was the last step)
1519 * we don't have to do it again, but otherwise print the final values.
1523 /* Write final value to log since we didn't do anything the last step */
1524 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1526 if (!do_ene || !do_log)
1528 /* Write final energy file entries */
1529 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1530 !do_log ? fplog : NULL, step, step, eprNORMAL,
1531 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1535 /* Print some stuff... */
1538 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1542 * For accurate normal mode calculation it is imperative that we
1543 * store the last conformation into the full precision binary trajectory.
1545 * However, we should only do it if we did NOT already write this step
1546 * above (which we did if do_x or do_f was true).
1548 do_x = !do_per_step(step, inputrec->nstxout);
1549 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1551 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1552 top_global, inputrec, step,
1553 s_min, state_global, f_global);
1555 fnormn = s_min->fnorm/sqrt(state_global->natoms);
1559 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1560 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1561 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1562 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1564 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1567 finish_em(cr, outf, walltime_accounting, wcycle);
1569 /* To print the actual number of steps we needed somewhere */
1570 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1573 } /* That's all folks */
1576 double do_lbfgs(FILE *fplog, t_commrec *cr,
1577 int nfile, const t_filenm fnm[],
1578 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
1579 int gmx_unused nstglobalcomm,
1580 gmx_vsite_t *vsite, gmx_constr_t constr,
1581 int gmx_unused stepout,
1582 t_inputrec *inputrec,
1583 gmx_mtop_t *top_global, t_fcdata *fcd,
1586 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1587 gmx_edsam_t gmx_unused ed,
1589 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1590 gmx_membed_t gmx_unused membed,
1591 real gmx_unused cpt_period, real gmx_unused max_hours,
1592 const char gmx_unused *deviceOptions,
1594 unsigned long gmx_unused Flags,
1595 gmx_walltime_accounting_t walltime_accounting)
1597 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1599 gmx_localtop_t *top;
1600 gmx_enerdata_t *enerd;
1602 gmx_global_stat_t gstat;
1605 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1606 double stepsize, gpa, gpb, gpc, tmp, minstep;
1607 real *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1608 real *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1609 real a, b, c, maxdelta, delta;
1610 real diag, Epot0, Epot, EpotA, EpotB, EpotC;
1611 real dgdx, dgdg, sq, yr, beta;
1613 gmx_bool converged, first;
1616 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1618 int start, end, number_steps;
1620 int i, k, m, n, nfmax, gf, step;
1627 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1632 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).");
1635 n = 3*state->natoms;
1636 nmaxcorr = inputrec->nbfgscorr;
1638 /* Allocate memory */
1639 /* Use pointers to real so we dont have to loop over both atoms and
1640 * dimensions all the time...
1641 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1642 * that point to the same memory.
1655 snew(rho, nmaxcorr);
1656 snew(alpha, nmaxcorr);
1659 for (i = 0; i < nmaxcorr; i++)
1665 for (i = 0; i < nmaxcorr; i++)
1674 init_em(fplog, LBFGS, cr, inputrec,
1675 state, top_global, &ems, &top, &f, &f_global,
1676 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1677 nfile, fnm, &outf, &mdebin, imdport, Flags);
1678 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1679 * so we free some memory again.
1684 xx = (real *)state->x;
1688 end = mdatoms->homenr;
1690 /* Print to log file */
1691 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1693 do_log = do_ene = do_x = do_f = TRUE;
1695 /* Max number of steps */
1696 number_steps = inputrec->nsteps;
1698 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1700 for (i = start; i < end; i++)
1702 if (mdatoms->cFREEZE)
1704 gf = mdatoms->cFREEZE[i];
1706 for (m = 0; m < DIM; m++)
1708 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1713 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1717 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1722 construct_vsites(vsite, state->x, 1, NULL,
1723 top->idef.iparams, top->idef.il,
1724 fr->ePBC, fr->bMolPBC, cr, state->box);
1727 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1728 /* do_force always puts the charge groups in the box and shifts again
1729 * We do not unshift, so molecules are always whole
1734 evaluate_energy(fplog, cr,
1735 top_global, &ems, top,
1736 inputrec, nrnb, wcycle, gstat,
1737 vsite, constr, fcd, graph, mdatoms, fr,
1738 mu_tot, enerd, vir, pres, -1, TRUE);
1743 /* Copy stuff to the energy bin for easy printing etc. */
1744 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1745 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
1746 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1748 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
1749 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1750 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1754 /* This is the starting energy */
1755 Epot = enerd->term[F_EPOT];
1761 /* Set the initial step.
1762 * since it will be multiplied by the non-normalized search direction
1763 * vector (force vector the first time), we scale it by the
1764 * norm of the force.
1769 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1770 fprintf(stderr, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1771 fprintf(stderr, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1772 fprintf(stderr, "\n");
1773 /* and copy to the log file too... */
1774 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1775 fprintf(fplog, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1776 fprintf(fplog, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1777 fprintf(fplog, "\n");
1781 for (i = 0; i < n; i++)
1785 dx[point][i] = ff[i]; /* Initial search direction */
1793 stepsize = 1.0/fnorm;
1795 /* Start the loop over BFGS steps.
1796 * Each successful step is counted, and we continue until
1797 * we either converge or reach the max number of steps.
1802 /* Set the gradient from the force */
1804 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1807 /* Write coordinates if necessary */
1808 do_x = do_per_step(step, inputrec->nstxout);
1809 do_f = do_per_step(step, inputrec->nstfout);
1814 mdof_flags |= MDOF_X;
1819 mdof_flags |= MDOF_F;
1824 mdof_flags |= MDOF_IMD;
1827 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1828 top_global, step, (real)step, state, state, f, f);
1830 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1832 /* pointer to current direction - point=0 first time here */
1835 /* calculate line gradient */
1836 for (gpa = 0, i = 0; i < n; i++)
1841 /* Calculate minimum allowed stepsize, before the average (norm)
1842 * relative change in coordinate is smaller than precision
1844 for (minstep = 0, i = 0; i < n; i++)
1854 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1856 if (stepsize < minstep)
1862 /* Store old forces and coordinates */
1863 for (i = 0; i < n; i++)
1872 for (i = 0; i < n; i++)
1877 /* Take a step downhill.
1878 * In theory, we should minimize the function along this direction.
1879 * That is quite possible, but it turns out to take 5-10 function evaluations
1880 * for each line. However, we dont really need to find the exact minimum -
1881 * it is much better to start a new BFGS step in a modified direction as soon
1882 * as we are close to it. This will save a lot of energy evaluations.
1884 * In practice, we just try to take a single step.
1885 * If it worked (i.e. lowered the energy), we increase the stepsize but
1886 * the continue straight to the next BFGS step without trying to find any minimum.
1887 * If it didn't work (higher energy), there must be a minimum somewhere between
1888 * the old position and the new one.
1890 * Due to the finite numerical accuracy, it turns out that it is a good idea
1891 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1892 * This leads to lower final energies in the tests I've done. / Erik
1897 c = a + stepsize; /* reference position along line is zero */
1899 /* Check stepsize first. We do not allow displacements
1900 * larger than emstep.
1906 for (i = 0; i < n; i++)
1909 if (delta > maxdelta)
1914 if (maxdelta > inputrec->em_stepsize)
1919 while (maxdelta > inputrec->em_stepsize);
1921 /* Take a trial step */
1922 for (i = 0; i < n; i++)
1924 xc[i] = lastx[i] + c*s[i];
1928 /* Calculate energy for the trial step */
1929 ems.s.x = (rvec *)xc;
1931 evaluate_energy(fplog, cr,
1932 top_global, &ems, top,
1933 inputrec, nrnb, wcycle, gstat,
1934 vsite, constr, fcd, graph, mdatoms, fr,
1935 mu_tot, enerd, vir, pres, step, FALSE);
1938 /* Calc derivative along line */
1939 for (gpc = 0, i = 0; i < n; i++)
1941 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1943 /* Sum the gradient along the line across CPUs */
1946 gmx_sumd(1, &gpc, cr);
1949 /* This is the max amount of increase in energy we tolerate */
1950 tmp = sqrt(GMX_REAL_EPS)*fabs(EpotA);
1952 /* Accept the step if the energy is lower, or if it is not significantly higher
1953 * and the line derivative is still negative.
1955 if (EpotC < EpotA || (gpc < 0 && EpotC < (EpotA+tmp)))
1958 /* Great, we found a better energy. Increase step for next iteration
1959 * if we are still going down, decrease it otherwise
1963 stepsize *= 1.618034; /* The golden section */
1967 stepsize *= 0.618034; /* 1/golden section */
1972 /* New energy is the same or higher. We will have to do some work
1973 * to find a smaller value in the interval. Take smaller step next time!
1976 stepsize *= 0.618034;
1979 /* OK, if we didn't find a lower value we will have to locate one now - there must
1980 * be one in the interval [a=0,c].
1981 * The same thing is valid here, though: Don't spend dozens of iterations to find
1982 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1983 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1985 * I also have a safeguard for potentially really patological functions so we never
1986 * take more than 20 steps before we give up ...
1988 * If we already found a lower value we just skip this step and continue to the update.
1997 /* Select a new trial point.
1998 * If the derivatives at points a & c have different sign we interpolate to zero,
1999 * otherwise just do a bisection.
2002 if (gpa < 0 && gpc > 0)
2004 b = a + gpa*(a-c)/(gpc-gpa);
2011 /* safeguard if interpolation close to machine accuracy causes errors:
2012 * never go outside the interval
2014 if (b <= a || b >= c)
2019 /* Take a trial step */
2020 for (i = 0; i < n; i++)
2022 xb[i] = lastx[i] + b*s[i];
2026 /* Calculate energy for the trial step */
2027 ems.s.x = (rvec *)xb;
2029 evaluate_energy(fplog, cr,
2030 top_global, &ems, top,
2031 inputrec, nrnb, wcycle, gstat,
2032 vsite, constr, fcd, graph, mdatoms, fr,
2033 mu_tot, enerd, vir, pres, step, FALSE);
2038 for (gpb = 0, i = 0; i < n; i++)
2040 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2043 /* Sum the gradient along the line across CPUs */
2046 gmx_sumd(1, &gpb, cr);
2049 /* Keep one of the intervals based on the value of the derivative at the new point */
2052 /* Replace c endpoint with b */
2056 /* swap coord pointers b/c */
2066 /* Replace a endpoint with b */
2070 /* swap coord pointers a/b */
2080 * Stop search as soon as we find a value smaller than the endpoints,
2081 * or if the tolerance is below machine precision.
2082 * Never run more than 20 steps, no matter what.
2086 while ((EpotB > EpotA || EpotB > EpotC) && (nminstep < 20));
2088 if (fabs(EpotB-Epot0) < GMX_REAL_EPS || nminstep >= 20)
2090 /* OK. We couldn't find a significantly lower energy.
2091 * If ncorr==0 this was steepest descent, and then we give up.
2092 * If not, reset memory to restart as steepest descent before quitting.
2104 /* Search in gradient direction */
2105 for (i = 0; i < n; i++)
2107 dx[point][i] = ff[i];
2109 /* Reset stepsize */
2110 stepsize = 1.0/fnorm;
2115 /* Select min energy state of A & C, put the best in xx/ff/Epot
2121 for (i = 0; i < n; i++)
2132 for (i = 0; i < n; i++)
2146 for (i = 0; i < n; i++)
2154 /* Update the memory information, and calculate a new
2155 * approximation of the inverse hessian
2158 /* Have new data in Epot, xx, ff */
2159 if (ncorr < nmaxcorr)
2164 for (i = 0; i < n; i++)
2166 dg[point][i] = lastf[i]-ff[i];
2167 dx[point][i] *= stepsize;
2172 for (i = 0; i < n; i++)
2174 dgdg += dg[point][i]*dg[point][i];
2175 dgdx += dg[point][i]*dx[point][i];
2180 rho[point] = 1.0/dgdx;
2183 if (point >= nmaxcorr)
2189 for (i = 0; i < n; i++)
2196 /* Recursive update. First go back over the memory points */
2197 for (k = 0; k < ncorr; k++)
2206 for (i = 0; i < n; i++)
2208 sq += dx[cp][i]*p[i];
2211 alpha[cp] = rho[cp]*sq;
2213 for (i = 0; i < n; i++)
2215 p[i] -= alpha[cp]*dg[cp][i];
2219 for (i = 0; i < n; i++)
2224 /* And then go forward again */
2225 for (k = 0; k < ncorr; k++)
2228 for (i = 0; i < n; i++)
2230 yr += p[i]*dg[cp][i];
2234 beta = alpha[cp]-beta;
2236 for (i = 0; i < n; i++)
2238 p[i] += beta*dx[cp][i];
2248 for (i = 0; i < n; i++)
2252 dx[point][i] = p[i];
2262 /* Test whether the convergence criterion is met */
2263 get_f_norm_max(cr, &(inputrec->opts), mdatoms, f, &fnorm, &fmax, &nfmax);
2265 /* Print it if necessary */
2270 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2271 step, Epot, fnorm/sqrt(state->natoms), fmax, nfmax+1);
2273 /* Store the new (lower) energies */
2274 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2275 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
2276 NULL, NULL, vir, pres, NULL, mu_tot, constr);
2277 do_log = do_per_step(step, inputrec->nstlog);
2278 do_ene = do_per_step(step, inputrec->nstenergy);
2281 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2283 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2284 do_log ? fplog : NULL, step, step, eprNORMAL,
2285 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2288 /* Send x and E to IMD client, if bIMD is TRUE. */
2289 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state->box, state->x, inputrec, 0, wcycle) && MASTER(cr))
2291 IMD_send_positions(inputrec->imd);
2294 /* Stop when the maximum force lies below tolerance.
2295 * If we have reached machine precision, converged is already set to true.
2298 converged = converged || (fmax < inputrec->em_tol);
2300 } /* End of the loop */
2302 /* IMD cleanup, if bIMD is TRUE. */
2303 IMD_finalize(inputrec->bIMD, inputrec->imd);
2307 step--; /* we never took that last step in this case */
2310 if (fmax > inputrec->em_tol)
2314 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2315 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2320 /* If we printed energy and/or logfile last step (which was the last step)
2321 * we don't have to do it again, but otherwise print the final values.
2323 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2325 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2327 if (!do_ene || !do_log) /* Write final energy file entries */
2329 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2330 !do_log ? fplog : NULL, step, step, eprNORMAL,
2331 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2334 /* Print some stuff... */
2337 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2341 * For accurate normal mode calculation it is imperative that we
2342 * store the last conformation into the full precision binary trajectory.
2344 * However, we should only do it if we did NOT already write this step
2345 * above (which we did if do_x or do_f was true).
2347 do_x = !do_per_step(step, inputrec->nstxout);
2348 do_f = !do_per_step(step, inputrec->nstfout);
2349 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2350 top_global, inputrec, step,
2355 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2356 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2357 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2358 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2360 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2363 finish_em(cr, outf, walltime_accounting, wcycle);
2365 /* To print the actual number of steps we needed somewhere */
2366 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2369 } /* That's all folks */
2372 double do_steep(FILE *fplog, t_commrec *cr,
2373 int nfile, const t_filenm fnm[],
2374 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2375 int gmx_unused nstglobalcomm,
2376 gmx_vsite_t *vsite, gmx_constr_t constr,
2377 int gmx_unused stepout,
2378 t_inputrec *inputrec,
2379 gmx_mtop_t *top_global, t_fcdata *fcd,
2380 t_state *state_global,
2382 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2383 gmx_edsam_t gmx_unused ed,
2385 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2386 gmx_membed_t gmx_unused membed,
2387 real gmx_unused cpt_period, real gmx_unused max_hours,
2388 const char gmx_unused *deviceOptions,
2390 unsigned long gmx_unused Flags,
2391 gmx_walltime_accounting_t walltime_accounting)
2393 const char *SD = "Steepest Descents";
2394 em_state_t *s_min, *s_try;
2396 gmx_localtop_t *top;
2397 gmx_enerdata_t *enerd;
2399 gmx_global_stat_t gstat;
2401 real stepsize, constepsize;
2405 gmx_bool bDone, bAbort, do_x, do_f;
2410 int steps_accepted = 0;
2414 s_min = init_em_state();
2415 s_try = init_em_state();
2417 /* Init em and store the local state in s_try */
2418 init_em(fplog, SD, cr, inputrec,
2419 state_global, top_global, s_try, &top, &f, &f_global,
2420 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2421 nfile, fnm, &outf, &mdebin, imdport, Flags);
2423 /* Print to log file */
2424 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2426 /* Set variables for stepsize (in nm). This is the largest
2427 * step that we are going to make in any direction.
2429 ustep = inputrec->em_stepsize;
2432 /* Max number of steps */
2433 nsteps = inputrec->nsteps;
2437 /* Print to the screen */
2438 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2442 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2445 /**** HERE STARTS THE LOOP ****
2446 * count is the counter for the number of steps
2447 * bDone will be TRUE when the minimization has converged
2448 * bAbort will be TRUE when nsteps steps have been performed or when
2449 * the stepsize becomes smaller than is reasonable for machine precision
2454 while (!bDone && !bAbort)
2456 bAbort = (nsteps >= 0) && (count == nsteps);
2458 /* set new coordinates, except for first step */
2461 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2462 s_min, stepsize, s_min->f, s_try,
2463 constr, top, nrnb, wcycle, count);
2466 evaluate_energy(fplog, cr,
2467 top_global, s_try, top,
2468 inputrec, nrnb, wcycle, gstat,
2469 vsite, constr, fcd, graph, mdatoms, fr,
2470 mu_tot, enerd, vir, pres, count, count == 0);
2474 print_ebin_header(fplog, count, count, s_try->s.lambda[efptFEP]);
2479 s_min->epot = s_try->epot + 1;
2482 /* Print it if necessary */
2487 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2488 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2489 (s_try->epot < s_min->epot) ? '\n' : '\r');
2492 if (s_try->epot < s_min->epot)
2494 /* Store the new (lower) energies */
2495 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2496 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2497 s_try->s.box, NULL, NULL, vir, pres, NULL, mu_tot, constr);
2499 /* Prepare IMD energy record, if bIMD is TRUE. */
2500 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2502 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2503 do_per_step(steps_accepted, inputrec->nstdisreout),
2504 do_per_step(steps_accepted, inputrec->nstorireout),
2505 fplog, count, count, eprNORMAL, TRUE,
2506 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2511 /* Now if the new energy is smaller than the previous...
2512 * or if this is the first step!
2513 * or if we did random steps!
2516 if ( (count == 0) || (s_try->epot < s_min->epot) )
2520 /* Test whether the convergence criterion is met... */
2521 bDone = (s_try->fmax < inputrec->em_tol);
2523 /* Copy the arrays for force, positions and energy */
2524 /* The 'Min' array always holds the coords and forces of the minimal
2526 swap_em_state(s_min, s_try);
2532 /* Write to trn, if necessary */
2533 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2534 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2535 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
2536 top_global, inputrec, count,
2537 s_min, state_global, f_global);
2541 /* If energy is not smaller make the step smaller... */
2544 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2546 /* Reload the old state */
2547 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2548 s_min, top, mdatoms, fr, vsite, constr,
2553 /* Determine new step */
2554 stepsize = ustep/s_min->fmax;
2556 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2558 if (count == nsteps || ustep < 1e-12)
2560 if (count == nsteps || ustep < 1e-6)
2565 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2566 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2571 /* Send IMD energies and positions, if bIMD is TRUE. */
2572 if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
2574 IMD_send_positions(inputrec->imd);
2578 } /* End of the loop */
2580 /* IMD cleanup, if bIMD is TRUE. */
2581 IMD_finalize(inputrec->bIMD, inputrec->imd);
2583 /* Print some data... */
2586 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2588 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2589 top_global, inputrec, count,
2590 s_min, state_global, f_global);
2592 fnormn = s_min->fnorm/sqrt(state_global->natoms);
2596 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2597 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2598 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2599 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2602 finish_em(cr, outf, walltime_accounting, wcycle);
2604 /* To print the actual number of steps we needed somewhere */
2605 inputrec->nsteps = count;
2607 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2610 } /* That's all folks */
2613 double do_nm(FILE *fplog, t_commrec *cr,
2614 int nfile, const t_filenm fnm[],
2615 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2616 int gmx_unused nstglobalcomm,
2617 gmx_vsite_t *vsite, gmx_constr_t constr,
2618 int gmx_unused stepout,
2619 t_inputrec *inputrec,
2620 gmx_mtop_t *top_global, t_fcdata *fcd,
2621 t_state *state_global,
2623 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2624 gmx_edsam_t gmx_unused ed,
2626 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2627 gmx_membed_t gmx_unused membed,
2628 real gmx_unused cpt_period, real gmx_unused max_hours,
2629 const char gmx_unused *deviceOptions,
2631 unsigned long gmx_unused Flags,
2632 gmx_walltime_accounting_t walltime_accounting)
2634 const char *NM = "Normal Mode Analysis";
2636 int natoms, atom, d;
2639 gmx_localtop_t *top;
2640 gmx_enerdata_t *enerd;
2642 gmx_global_stat_t gstat;
2644 real t, t0, lambda, lam0;
2649 gmx_bool bSparse; /* use sparse matrix storage format */
2651 gmx_sparsematrix_t * sparse_matrix = NULL;
2652 real * full_matrix = NULL;
2653 em_state_t * state_work;
2655 /* added with respect to mdrun */
2656 int i, j, k, row, col;
2657 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2663 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2666 state_work = init_em_state();
2668 /* Init em and store the local state in state_minimum */
2669 init_em(fplog, NM, cr, inputrec,
2670 state_global, top_global, state_work, &top,
2672 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2673 nfile, fnm, &outf, NULL, imdport, Flags);
2675 natoms = top_global->natoms;
2683 "NOTE: This version of Gromacs has been compiled in single precision,\n"
2684 " which MIGHT not be accurate enough for normal mode analysis.\n"
2685 " Gromacs now uses sparse matrix storage, so the memory requirements\n"
2686 " are fairly modest even if you recompile in double precision.\n\n");
2690 /* Check if we can/should use sparse storage format.
2692 * Sparse format is only useful when the Hessian itself is sparse, which it
2693 * will be when we use a cutoff.
2694 * For small systems (n<1000) it is easier to always use full matrix format, though.
2696 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2698 md_print_info(cr, fplog, "Non-cutoff electrostatics used, forcing full Hessian format.\n");
2701 else if (top_global->natoms < 1000)
2703 md_print_info(cr, fplog, "Small system size (N=%d), using full Hessian format.\n", top_global->natoms);
2708 md_print_info(cr, fplog, "Using compressed symmetric sparse Hessian format.\n");
2714 sz = DIM*top_global->natoms;
2716 fprintf(stderr, "Allocating Hessian memory...\n\n");
2720 sparse_matrix = gmx_sparsematrix_init(sz);
2721 sparse_matrix->compressed_symmetric = TRUE;
2725 snew(full_matrix, sz*sz);
2729 /* Initial values */
2730 t0 = inputrec->init_t;
2731 lam0 = inputrec->fepvals->init_lambda;
2739 /* Write start time and temperature */
2740 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2742 /* fudge nr of steps to nr of atoms */
2743 inputrec->nsteps = natoms*2;
2747 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2748 *(top_global->name), (int)inputrec->nsteps);
2751 nnodes = cr->nnodes;
2753 /* Make evaluate_energy do a single node force calculation */
2755 evaluate_energy(fplog, cr,
2756 top_global, state_work, top,
2757 inputrec, nrnb, wcycle, gstat,
2758 vsite, constr, fcd, graph, mdatoms, fr,
2759 mu_tot, enerd, vir, pres, -1, TRUE);
2760 cr->nnodes = nnodes;
2762 /* if forces are not small, warn user */
2763 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, state_work);
2765 md_print_info(cr, fplog, "Maximum force:%12.5e\n", state_work->fmax);
2766 if (state_work->fmax > 1.0e-3)
2768 md_print_info(cr, fplog,
2769 "The force is probably not small enough to "
2770 "ensure that you are at a minimum.\n"
2771 "Be aware that negative eigenvalues may occur\n"
2772 "when the resulting matrix is diagonalized.\n\n");
2775 /***********************************************************
2777 * Loop over all pairs in matrix
2779 * do_force called twice. Once with positive and
2780 * once with negative displacement
2782 ************************************************************/
2784 /* Steps are divided one by one over the nodes */
2785 for (atom = cr->nodeid; atom < natoms; atom += nnodes)
2788 for (d = 0; d < DIM; d++)
2790 x_min = state_work->s.x[atom][d];
2792 state_work->s.x[atom][d] = x_min - der_range;
2794 /* Make evaluate_energy do a single node force calculation */
2796 evaluate_energy(fplog, cr,
2797 top_global, state_work, top,
2798 inputrec, nrnb, wcycle, gstat,
2799 vsite, constr, fcd, graph, mdatoms, fr,
2800 mu_tot, enerd, vir, pres, atom*2, FALSE);
2802 for (i = 0; i < natoms; i++)
2804 copy_rvec(state_work->f[i], fneg[i]);
2807 state_work->s.x[atom][d] = x_min + der_range;
2809 evaluate_energy(fplog, cr,
2810 top_global, state_work, top,
2811 inputrec, nrnb, wcycle, gstat,
2812 vsite, constr, fcd, graph, mdatoms, fr,
2813 mu_tot, enerd, vir, pres, atom*2+1, FALSE);
2814 cr->nnodes = nnodes;
2816 /* x is restored to original */
2817 state_work->s.x[atom][d] = x_min;
2819 for (j = 0; j < natoms; j++)
2821 for (k = 0; (k < DIM); k++)
2824 -(state_work->f[j][k] - fneg[j][k])/(2*der_range);
2832 #define mpi_type MPI_DOUBLE
2834 #define mpi_type MPI_FLOAT
2836 MPI_Send(dfdx[0], natoms*DIM, mpi_type, MASTERNODE(cr), cr->nodeid,
2837 cr->mpi_comm_mygroup);
2842 for (node = 0; (node < nnodes && atom+node < natoms); node++)
2848 MPI_Recv(dfdx[0], natoms*DIM, mpi_type, node, node,
2849 cr->mpi_comm_mygroup, &stat);
2854 row = (atom + node)*DIM + d;
2856 for (j = 0; j < natoms; j++)
2858 for (k = 0; k < DIM; k++)
2864 if (col >= row && dfdx[j][k] != 0.0)
2866 gmx_sparsematrix_increment_value(sparse_matrix,
2867 row, col, dfdx[j][k]);
2872 full_matrix[row*sz+col] = dfdx[j][k];
2879 if (bVerbose && fplog)
2884 /* write progress */
2885 if (MASTER(cr) && bVerbose)
2887 fprintf(stderr, "\rFinished step %d out of %d",
2888 min(atom+nnodes, natoms), natoms);
2895 fprintf(stderr, "\n\nWriting Hessian...\n");
2896 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2899 finish_em(cr, outf, walltime_accounting, wcycle);
2901 walltime_accounting_set_nsteps_done(walltime_accounting, natoms*2);