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.
52 #include "gromacs/random/random.h"
54 #include "gmx_fatal.h"
65 #include "md_support.h"
71 #include "mtop_util.h"
74 #include "gmx_omp_nthreads.h"
75 #include "md_logging.h"
77 #include "gromacs/fileio/confio.h"
78 #include "gromacs/fileio/trajectory_writing.h"
79 #include "gromacs/linearalgebra/mtxio.h"
80 #include "gromacs/linearalgebra/sparsematrix.h"
81 #include "gromacs/timing/wallcycle.h"
82 #include "gromacs/timing/walltime_accounting.h"
93 static em_state_t *init_em_state()
99 /* does this need to be here? Should the array be declared differently (staticaly)in the state definition? */
100 snew(ems->s.lambda, efptNR);
105 static void print_em_start(FILE *fplog,
107 gmx_walltime_accounting_t walltime_accounting,
108 gmx_wallcycle_t wcycle,
111 walltime_accounting_start(walltime_accounting);
112 wallcycle_start(wcycle, ewcRUN);
113 print_start(fplog, cr, walltime_accounting, name);
115 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
116 gmx_wallcycle_t wcycle)
118 wallcycle_stop(wcycle, ewcRUN);
120 walltime_accounting_end(walltime_accounting);
123 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
126 fprintf(out, "%s:\n", minimizer);
127 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
128 fprintf(out, " Number of steps = %12d\n", nsteps);
131 static void warn_step(FILE *fp, real ftol, gmx_bool bLastStep, gmx_bool bConstrain)
137 "\nEnergy minimization reached the maximum number "
138 "of steps before the forces reached the requested "
139 "precision Fmax < %g.\n", ftol);
144 "\nEnergy minimization has stopped, but the forces have "
145 "not converged to the requested precision Fmax < %g (which "
146 "may not be possible for your system). It stopped "
147 "because the algorithm tried to make a new step whose size "
148 "was too small, or there was no change in the energy since "
149 "last step. Either way, we regard the minimization as "
150 "converged to within the available machine precision, "
151 "given your starting configuration and EM parameters.\n%s%s",
153 sizeof(real) < sizeof(double) ?
154 "\nDouble precision normally gives you higher accuracy, but "
155 "this is often not needed for preparing to run molecular "
159 "You might need to increase your constraint accuracy, or turn\n"
160 "off constraints altogether (set constraints = none in mdp file)\n" :
163 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
168 static void print_converged(FILE *fp, const char *alg, real ftol,
169 gmx_int64_t count, gmx_bool bDone, gmx_int64_t nsteps,
170 real epot, real fmax, int nfmax, real fnorm)
172 char buf[STEPSTRSIZE];
176 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
177 alg, ftol, gmx_step_str(count, buf));
179 else if (count < nsteps)
181 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
182 "but did not reach the requested Fmax < %g.\n",
183 alg, gmx_step_str(count, buf), ftol);
187 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
188 alg, ftol, gmx_step_str(count, buf));
192 fprintf(fp, "Potential Energy = %21.14e\n", epot);
193 fprintf(fp, "Maximum force = %21.14e on atom %d\n", fmax, nfmax+1);
194 fprintf(fp, "Norm of force = %21.14e\n", fnorm);
196 fprintf(fp, "Potential Energy = %14.7e\n", epot);
197 fprintf(fp, "Maximum force = %14.7e on atom %d\n", fmax, nfmax+1);
198 fprintf(fp, "Norm of force = %14.7e\n", fnorm);
202 static void get_f_norm_max(t_commrec *cr,
203 t_grpopts *opts, t_mdatoms *mdatoms, rvec *f,
204 real *fnorm, real *fmax, int *a_fmax)
207 real fmax2, fmax2_0, fam;
208 int la_max, a_max, start, end, i, m, gf;
210 /* This routine finds the largest force and returns it.
211 * On parallel machines the global max is taken.
217 start = mdatoms->start;
218 end = mdatoms->homenr + start;
219 if (mdatoms->cFREEZE)
221 for (i = start; i < end; i++)
223 gf = mdatoms->cFREEZE[i];
225 for (m = 0; m < DIM; m++)
227 if (!opts->nFreeze[gf][m])
242 for (i = start; i < end; i++)
254 if (la_max >= 0 && DOMAINDECOMP(cr))
256 a_max = cr->dd->gatindex[la_max];
264 snew(sum, 2*cr->nnodes+1);
265 sum[2*cr->nodeid] = fmax2;
266 sum[2*cr->nodeid+1] = a_max;
267 sum[2*cr->nnodes] = fnorm2;
268 gmx_sumd(2*cr->nnodes+1, sum, cr);
269 fnorm2 = sum[2*cr->nnodes];
270 /* Determine the global maximum */
271 for (i = 0; i < cr->nnodes; i++)
273 if (sum[2*i] > fmax2)
276 a_max = (int)(sum[2*i+1] + 0.5);
284 *fnorm = sqrt(fnorm2);
296 static void get_state_f_norm_max(t_commrec *cr,
297 t_grpopts *opts, t_mdatoms *mdatoms,
300 get_f_norm_max(cr, opts, mdatoms, ems->f, &ems->fnorm, &ems->fmax, &ems->a_fmax);
303 void init_em(FILE *fplog, const char *title,
304 t_commrec *cr, t_inputrec *ir,
305 t_state *state_global, gmx_mtop_t *top_global,
306 em_state_t *ems, gmx_localtop_t **top,
307 rvec **f, rvec **f_global,
308 t_nrnb *nrnb, rvec mu_tot,
309 t_forcerec *fr, gmx_enerdata_t **enerd,
310 t_graph **graph, t_mdatoms *mdatoms, gmx_global_stat_t *gstat,
311 gmx_vsite_t *vsite, gmx_constr_t constr,
312 int nfile, const t_filenm fnm[],
313 gmx_mdoutf_t *outf, t_mdebin **mdebin)
315 int start, homenr, i;
320 fprintf(fplog, "Initiating %s\n", title);
323 state_global->ngtc = 0;
325 /* Initialize lambda variables */
326 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, NULL);
330 if (DOMAINDECOMP(cr))
332 *top = dd_init_local_top(top_global);
334 dd_init_local_state(cr->dd, state_global, &ems->s);
338 /* Distribute the charge groups over the nodes from the master node */
339 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
340 state_global, top_global, ir,
341 &ems->s, &ems->f, mdatoms, *top,
342 fr, vsite, NULL, constr,
344 dd_store_state(cr->dd, &ems->s);
348 snew(*f_global, top_global->natoms);
358 snew(*f, top_global->natoms);
360 /* Just copy the state */
361 ems->s = *state_global;
362 snew(ems->s.x, ems->s.nalloc);
363 snew(ems->f, ems->s.nalloc);
364 for (i = 0; i < state_global->natoms; i++)
366 copy_rvec(state_global->x[i], ems->s.x[i]);
368 copy_mat(state_global->box, ems->s.box);
370 if (PAR(cr) && ir->eI != eiNM)
372 /* Initialize the particle decomposition and split the topology */
373 *top = split_system(fplog, top_global, ir, cr);
375 pd_cg_range(cr, &fr->cg0, &fr->hcg);
379 *top = gmx_mtop_generate_local_top(top_global, ir);
383 forcerec_set_excl_load(fr, *top, cr);
385 setup_bonded_threading(fr, &(*top)->idef);
387 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
389 *graph = mk_graph(fplog, &((*top)->idef), 0, top_global->natoms, FALSE, FALSE);
398 pd_at_range(cr, &start, &homenr);
404 homenr = top_global->natoms;
406 atoms2md(top_global, ir, 0, NULL, start, homenr, mdatoms);
407 update_mdatoms(mdatoms, state_global->lambda[efptFEP]);
411 set_vsite_top(vsite, *top, mdatoms, cr);
417 if (ir->eConstrAlg == econtSHAKE &&
418 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
420 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
421 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
424 if (!DOMAINDECOMP(cr))
426 set_constraints(constr, *top, ir, mdatoms, cr);
429 if (!ir->bContinuation)
431 /* Constrain the starting coordinates */
433 constrain(PAR(cr) ? NULL : fplog, TRUE, TRUE, constr, &(*top)->idef,
434 ir, NULL, cr, -1, 0, mdatoms,
435 ems->s.x, ems->s.x, NULL, fr->bMolPBC, ems->s.box,
436 ems->s.lambda[efptFEP], &dvdl_constr,
437 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
443 *gstat = global_stat_init(ir);
446 *outf = init_mdoutf(nfile, fnm, 0, cr, ir, top_global, NULL);
449 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
454 /* Init bin for energy stuff */
455 *mdebin = init_mdebin(mdoutf_get_fp_ene(*outf), top_global, ir, NULL);
459 calc_shifts(ems->s.box, fr->shift_vec);
462 static void finish_em(t_commrec *cr, gmx_mdoutf_t outf,
463 gmx_walltime_accounting_t walltime_accounting,
464 gmx_wallcycle_t wcycle)
466 if (!(cr->duty & DUTY_PME))
468 /* Tell the PME only node to finish */
469 gmx_pme_send_finish(cr);
474 em_time_end(walltime_accounting, wcycle);
477 static void swap_em_state(em_state_t *ems1, em_state_t *ems2)
486 static void copy_em_coords(em_state_t *ems, t_state *state)
490 for (i = 0; (i < state->natoms); i++)
492 copy_rvec(ems->s.x[i], state->x[i]);
496 static void write_em_traj(FILE *fplog, t_commrec *cr,
498 gmx_bool bX, gmx_bool bF, const char *confout,
499 gmx_mtop_t *top_global,
500 t_inputrec *ir, gmx_int64_t step,
502 t_state *state_global, rvec *f_global)
506 if ((bX || bF || confout != NULL) && !DOMAINDECOMP(cr))
508 copy_em_coords(state, state_global);
515 mdof_flags |= MDOF_X;
519 mdof_flags |= MDOF_F;
521 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
522 top_global, step, (double)step,
523 &state->s, state_global, state->f, f_global);
525 if (confout != NULL && MASTER(cr))
527 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
529 /* Make molecules whole only for confout writing */
530 do_pbc_mtop(fplog, ir->ePBC, state_global->box, top_global,
534 write_sto_conf_mtop(confout,
535 *top_global->name, top_global,
536 state_global->x, NULL, ir->ePBC, state_global->box);
540 static void do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
542 em_state_t *ems1, real a, rvec *f, em_state_t *ems2,
543 gmx_constr_t constr, gmx_localtop_t *top,
544 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
557 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
559 gmx_incons("state mismatch in do_em_step");
562 s2->flags = s1->flags;
564 if (s2->nalloc != s1->nalloc)
566 s2->nalloc = s1->nalloc;
567 srenew(s2->x, s1->nalloc);
568 srenew(ems2->f, s1->nalloc);
569 if (s2->flags & (1<<estCGP))
571 srenew(s2->cg_p, s1->nalloc);
575 s2->natoms = s1->natoms;
576 copy_mat(s1->box, s2->box);
577 /* Copy free energy state */
578 for (i = 0; i < efptNR; i++)
580 s2->lambda[i] = s1->lambda[i];
582 copy_mat(s1->box, s2->box);
585 end = md->start + md->homenr;
590 #pragma omp parallel num_threads(gmx_omp_nthreads_get(emntUpdate))
595 #pragma omp for schedule(static) nowait
596 for (i = start; i < end; i++)
602 for (m = 0; m < DIM; m++)
604 if (ir->opts.nFreeze[gf][m])
610 x2[i][m] = x1[i][m] + a*f[i][m];
615 if (s2->flags & (1<<estCGP))
617 /* Copy the CG p vector */
620 #pragma omp for schedule(static) nowait
621 for (i = start; i < end; i++)
623 copy_rvec(x1[i], x2[i]);
627 if (DOMAINDECOMP(cr))
629 s2->ddp_count = s1->ddp_count;
630 if (s2->cg_gl_nalloc < s1->cg_gl_nalloc)
633 s2->cg_gl_nalloc = s1->cg_gl_nalloc;
634 srenew(s2->cg_gl, s2->cg_gl_nalloc);
637 s2->ncg_gl = s1->ncg_gl;
638 #pragma omp for schedule(static) nowait
639 for (i = 0; i < s2->ncg_gl; i++)
641 s2->cg_gl[i] = s1->cg_gl[i];
643 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
649 wallcycle_start(wcycle, ewcCONSTR);
651 constrain(NULL, TRUE, TRUE, constr, &top->idef,
652 ir, NULL, cr, count, 0, md,
653 s1->x, s2->x, NULL, bMolPBC, s2->box,
654 s2->lambda[efptBONDED], &dvdl_constr,
655 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
656 wallcycle_stop(wcycle, ewcCONSTR);
660 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
661 gmx_mtop_t *top_global, t_inputrec *ir,
662 em_state_t *ems, gmx_localtop_t *top,
663 t_mdatoms *mdatoms, t_forcerec *fr,
664 gmx_vsite_t *vsite, gmx_constr_t constr,
665 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
667 /* Repartition the domain decomposition */
668 wallcycle_start(wcycle, ewcDOMDEC);
669 dd_partition_system(fplog, step, cr, FALSE, 1,
670 NULL, top_global, ir,
672 mdatoms, top, fr, vsite, NULL, constr,
673 nrnb, wcycle, FALSE);
674 dd_store_state(cr->dd, &ems->s);
675 wallcycle_stop(wcycle, ewcDOMDEC);
678 static void evaluate_energy(FILE *fplog, t_commrec *cr,
679 gmx_mtop_t *top_global,
680 em_state_t *ems, gmx_localtop_t *top,
681 t_inputrec *inputrec,
682 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
683 gmx_global_stat_t gstat,
684 gmx_vsite_t *vsite, gmx_constr_t constr,
686 t_graph *graph, t_mdatoms *mdatoms,
687 t_forcerec *fr, rvec mu_tot,
688 gmx_enerdata_t *enerd, tensor vir, tensor pres,
689 gmx_int64_t count, gmx_bool bFirst)
694 tensor force_vir, shake_vir, ekin;
695 real dvdl_constr, prescorr, enercorr, dvdlcorr;
698 /* Set the time to the initial time, the time does not change during EM */
699 t = inputrec->init_t;
702 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
704 /* This the first state or an old state used before the last ns */
710 if (inputrec->nstlist > 0)
714 else if (inputrec->nstlist == -1)
716 nabnsb = natoms_beyond_ns_buffer(inputrec, fr, &top->cgs, NULL, ems->s.x);
719 gmx_sumi(1, &nabnsb, cr);
727 construct_vsites(vsite, ems->s.x, 1, NULL,
728 top->idef.iparams, top->idef.il,
729 fr->ePBC, fr->bMolPBC, graph, cr, ems->s.box);
732 if (DOMAINDECOMP(cr))
736 /* Repartition the domain decomposition */
737 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
738 ems, top, mdatoms, fr, vsite, constr,
743 /* Calc force & energy on new trial position */
744 /* do_force always puts the charge groups in the box and shifts again
745 * We do not unshift, so molecules are always whole in congrad.c
747 do_force(fplog, cr, inputrec,
748 count, nrnb, wcycle, top, &top_global->groups,
749 ems->s.box, ems->s.x, &ems->s.hist,
750 ems->f, force_vir, mdatoms, enerd, fcd,
751 ems->s.lambda, graph, fr, vsite, mu_tot, t, NULL, NULL, TRUE,
752 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
753 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
754 (bNS ? GMX_FORCE_NS | GMX_FORCE_DO_LR : 0));
756 /* Clear the unused shake virial and pressure */
757 clear_mat(shake_vir);
760 /* Communicate stuff when parallel */
761 if (PAR(cr) && inputrec->eI != eiNM)
763 wallcycle_start(wcycle, ewcMoveE);
765 global_stat(fplog, gstat, cr, enerd, force_vir, shake_vir, mu_tot,
766 inputrec, NULL, NULL, NULL, 1, &terminate,
767 top_global, &ems->s, FALSE,
773 wallcycle_stop(wcycle, ewcMoveE);
776 /* Calculate long range corrections to pressure and energy */
777 calc_dispcorr(fplog, inputrec, fr, count, top_global->natoms, ems->s.box, ems->s.lambda[efptVDW],
778 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
779 enerd->term[F_DISPCORR] = enercorr;
780 enerd->term[F_EPOT] += enercorr;
781 enerd->term[F_PRES] += prescorr;
782 enerd->term[F_DVDL] += dvdlcorr;
784 ems->epot = enerd->term[F_EPOT];
788 /* Project out the constraint components of the force */
789 wallcycle_start(wcycle, ewcCONSTR);
791 constrain(NULL, FALSE, FALSE, constr, &top->idef,
792 inputrec, NULL, cr, count, 0, mdatoms,
793 ems->s.x, ems->f, ems->f, fr->bMolPBC, ems->s.box,
794 ems->s.lambda[efptBONDED], &dvdl_constr,
795 NULL, &shake_vir, nrnb, econqForceDispl, FALSE, 0, 0);
796 if (fr->bSepDVDL && fplog)
798 gmx_print_sepdvdl(fplog, "Constraints", t, dvdl_constr);
800 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
801 m_add(force_vir, shake_vir, vir);
802 wallcycle_stop(wcycle, ewcCONSTR);
806 copy_mat(force_vir, vir);
810 enerd->term[F_PRES] =
811 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
813 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
815 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
817 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
821 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
823 em_state_t *s_min, em_state_t *s_b)
827 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
829 unsigned char *grpnrFREEZE;
833 fprintf(debug, "Doing reorder_partsum\n");
839 cgs_gl = dd_charge_groups_global(cr->dd);
840 index = cgs_gl->index;
842 /* Collect fm in a global vector fmg.
843 * This conflicts with the spirit of domain decomposition,
844 * but to fully optimize this a much more complicated algorithm is required.
846 snew(fmg, mtop->natoms);
848 ncg = s_min->s.ncg_gl;
849 cg_gl = s_min->s.cg_gl;
851 for (c = 0; c < ncg; c++)
856 for (a = a0; a < a1; a++)
858 copy_rvec(fm[i], fmg[a]);
862 gmx_sum(mtop->natoms*3, fmg[0], cr);
864 /* Now we will determine the part of the sum for the cgs in state s_b */
866 cg_gl = s_b->s.cg_gl;
870 grpnrFREEZE = mtop->groups.grpnr[egcFREEZE];
871 for (c = 0; c < ncg; c++)
876 for (a = a0; a < a1; a++)
878 if (mdatoms->cFREEZE && grpnrFREEZE)
882 for (m = 0; m < DIM; m++)
884 if (!opts->nFreeze[gf][m])
886 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
898 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
900 em_state_t *s_min, em_state_t *s_b)
906 /* This is just the classical Polak-Ribiere calculation of beta;
907 * it looks a bit complicated since we take freeze groups into account,
908 * and might have to sum it in parallel runs.
911 if (!DOMAINDECOMP(cr) ||
912 (s_min->s.ddp_count == cr->dd->ddp_count &&
913 s_b->s.ddp_count == cr->dd->ddp_count))
919 /* This part of code can be incorrect with DD,
920 * since the atom ordering in s_b and s_min might differ.
922 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
924 if (mdatoms->cFREEZE)
926 gf = mdatoms->cFREEZE[i];
928 for (m = 0; m < DIM; m++)
930 if (!opts->nFreeze[gf][m])
932 sum += (fb[i][m] - fm[i][m])*fb[i][m];
939 /* We need to reorder cgs while summing */
940 sum = reorder_partsum(cr, opts, mdatoms, mtop, s_min, s_b);
944 gmx_sumd(1, &sum, cr);
947 return sum/sqr(s_min->fnorm);
950 double do_cg(FILE *fplog, t_commrec *cr,
951 int nfile, const t_filenm fnm[],
952 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
953 int gmx_unused nstglobalcomm,
954 gmx_vsite_t *vsite, gmx_constr_t constr,
955 int gmx_unused stepout,
956 t_inputrec *inputrec,
957 gmx_mtop_t *top_global, t_fcdata *fcd,
958 t_state *state_global,
960 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
961 gmx_edsam_t gmx_unused ed,
963 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
964 gmx_membed_t gmx_unused membed,
965 real gmx_unused cpt_period, real gmx_unused max_hours,
966 const char gmx_unused *deviceOptions,
967 unsigned long gmx_unused Flags,
968 gmx_walltime_accounting_t walltime_accounting)
970 const char *CG = "Polak-Ribiere Conjugate Gradients";
972 em_state_t *s_min, *s_a, *s_b, *s_c;
974 gmx_enerdata_t *enerd;
976 gmx_global_stat_t gstat;
978 rvec *f_global, *p, *sf, *sfm;
979 double gpa, gpb, gpc, tmp, sum[2], minstep;
982 real a, b, c, beta = 0.0;
986 gmx_bool converged, foundlower;
988 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
990 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
992 int i, m, gf, step, nminstep;
997 s_min = init_em_state();
998 s_a = init_em_state();
999 s_b = init_em_state();
1000 s_c = init_em_state();
1002 /* Init em and store the local state in s_min */
1003 init_em(fplog, CG, cr, inputrec,
1004 state_global, top_global, s_min, &top, &f, &f_global,
1005 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1006 nfile, fnm, &outf, &mdebin);
1008 /* Print to log file */
1009 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1011 /* Max number of steps */
1012 number_steps = inputrec->nsteps;
1016 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1020 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1023 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1024 /* do_force always puts the charge groups in the box and shifts again
1025 * We do not unshift, so molecules are always whole in congrad.c
1027 evaluate_energy(fplog, cr,
1028 top_global, s_min, top,
1029 inputrec, nrnb, wcycle, gstat,
1030 vsite, constr, fcd, graph, mdatoms, fr,
1031 mu_tot, enerd, vir, pres, -1, TRUE);
1036 /* Copy stuff to the energy bin for easy printing etc. */
1037 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1038 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1039 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1041 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1042 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1043 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1047 /* Estimate/guess the initial stepsize */
1048 stepsize = inputrec->em_stepsize/s_min->fnorm;
1052 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1053 s_min->fmax, s_min->a_fmax+1);
1054 fprintf(stderr, " F-Norm = %12.5e\n",
1055 s_min->fnorm/sqrt(state_global->natoms));
1056 fprintf(stderr, "\n");
1057 /* and copy to the log file too... */
1058 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1059 s_min->fmax, s_min->a_fmax+1);
1060 fprintf(fplog, " F-Norm = %12.5e\n",
1061 s_min->fnorm/sqrt(state_global->natoms));
1062 fprintf(fplog, "\n");
1064 /* Start the loop over CG steps.
1065 * Each successful step is counted, and we continue until
1066 * we either converge or reach the max number of steps.
1069 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1072 /* start taking steps in a new direction
1073 * First time we enter the routine, beta=0, and the direction is
1074 * simply the negative gradient.
1077 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1082 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1084 if (mdatoms->cFREEZE)
1086 gf = mdatoms->cFREEZE[i];
1088 for (m = 0; m < DIM; m++)
1090 if (!inputrec->opts.nFreeze[gf][m])
1092 p[i][m] = sf[i][m] + beta*p[i][m];
1093 gpa -= p[i][m]*sf[i][m];
1094 /* f is negative gradient, thus the sign */
1103 /* Sum the gradient along the line across CPUs */
1106 gmx_sumd(1, &gpa, cr);
1109 /* Calculate the norm of the search vector */
1110 get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1112 /* Just in case stepsize reaches zero due to numerical precision... */
1115 stepsize = inputrec->em_stepsize/pnorm;
1119 * Double check the value of the derivative in the search direction.
1120 * If it is positive it must be due to the old information in the
1121 * CG formula, so just remove that and start over with beta=0.
1122 * This corresponds to a steepest descent step.
1127 step--; /* Don't count this step since we are restarting */
1128 continue; /* Go back to the beginning of the big for-loop */
1131 /* Calculate minimum allowed stepsize, before the average (norm)
1132 * relative change in coordinate is smaller than precision
1135 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1137 for (m = 0; m < DIM; m++)
1139 tmp = fabs(s_min->s.x[i][m]);
1148 /* Add up from all CPUs */
1151 gmx_sumd(1, &minstep, cr);
1154 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1156 if (stepsize < minstep)
1162 /* Write coordinates if necessary */
1163 do_x = do_per_step(step, inputrec->nstxout);
1164 do_f = do_per_step(step, inputrec->nstfout);
1166 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1167 top_global, inputrec, step,
1168 s_min, state_global, f_global);
1170 /* Take a step downhill.
1171 * In theory, we should minimize the function along this direction.
1172 * That is quite possible, but it turns out to take 5-10 function evaluations
1173 * for each line. However, we dont really need to find the exact minimum -
1174 * it is much better to start a new CG step in a modified direction as soon
1175 * as we are close to it. This will save a lot of energy evaluations.
1177 * In practice, we just try to take a single step.
1178 * If it worked (i.e. lowered the energy), we increase the stepsize but
1179 * the continue straight to the next CG step without trying to find any minimum.
1180 * If it didn't work (higher energy), there must be a minimum somewhere between
1181 * the old position and the new one.
1183 * Due to the finite numerical accuracy, it turns out that it is a good idea
1184 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1185 * This leads to lower final energies in the tests I've done. / Erik
1187 s_a->epot = s_min->epot;
1189 c = a + stepsize; /* reference position along line is zero */
1191 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1193 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1194 s_min, top, mdatoms, fr, vsite, constr,
1198 /* Take a trial step (new coords in s_c) */
1199 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1200 constr, top, nrnb, wcycle, -1);
1203 /* Calculate energy for the trial step */
1204 evaluate_energy(fplog, cr,
1205 top_global, s_c, top,
1206 inputrec, nrnb, wcycle, gstat,
1207 vsite, constr, fcd, graph, mdatoms, fr,
1208 mu_tot, enerd, vir, pres, -1, FALSE);
1210 /* Calc derivative along line */
1214 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1216 for (m = 0; m < DIM; m++)
1218 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1221 /* Sum the gradient along the line across CPUs */
1224 gmx_sumd(1, &gpc, cr);
1227 /* This is the max amount of increase in energy we tolerate */
1228 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1230 /* Accept the step if the energy is lower, or if it is not significantly higher
1231 * and the line derivative is still negative.
1233 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1236 /* Great, we found a better energy. Increase step for next iteration
1237 * if we are still going down, decrease it otherwise
1241 stepsize *= 1.618034; /* The golden section */
1245 stepsize *= 0.618034; /* 1/golden section */
1250 /* New energy is the same or higher. We will have to do some work
1251 * to find a smaller value in the interval. Take smaller step next time!
1254 stepsize *= 0.618034;
1260 /* OK, if we didn't find a lower value we will have to locate one now - there must
1261 * be one in the interval [a=0,c].
1262 * The same thing is valid here, though: Don't spend dozens of iterations to find
1263 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1264 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1266 * I also have a safeguard for potentially really patological functions so we never
1267 * take more than 20 steps before we give up ...
1269 * If we already found a lower value we just skip this step and continue to the update.
1277 /* Select a new trial point.
1278 * If the derivatives at points a & c have different sign we interpolate to zero,
1279 * otherwise just do a bisection.
1281 if (gpa < 0 && gpc > 0)
1283 b = a + gpa*(a-c)/(gpc-gpa);
1290 /* safeguard if interpolation close to machine accuracy causes errors:
1291 * never go outside the interval
1293 if (b <= a || b >= c)
1298 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1300 /* Reload the old state */
1301 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1302 s_min, top, mdatoms, fr, vsite, constr,
1306 /* Take a trial step to this new point - new coords in s_b */
1307 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1308 constr, top, nrnb, wcycle, -1);
1311 /* Calculate energy for the trial step */
1312 evaluate_energy(fplog, cr,
1313 top_global, s_b, top,
1314 inputrec, nrnb, wcycle, gstat,
1315 vsite, constr, fcd, graph, mdatoms, fr,
1316 mu_tot, enerd, vir, pres, -1, FALSE);
1318 /* p does not change within a step, but since the domain decomposition
1319 * might change, we have to use cg_p of s_b here.
1324 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1326 for (m = 0; m < DIM; m++)
1328 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1331 /* Sum the gradient along the line across CPUs */
1334 gmx_sumd(1, &gpb, cr);
1339 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1340 s_a->epot, s_b->epot, s_c->epot, gpb);
1343 epot_repl = s_b->epot;
1345 /* Keep one of the intervals based on the value of the derivative at the new point */
1348 /* Replace c endpoint with b */
1349 swap_em_state(s_b, s_c);
1355 /* Replace a endpoint with b */
1356 swap_em_state(s_b, s_a);
1362 * Stop search as soon as we find a value smaller than the endpoints.
1363 * Never run more than 20 steps, no matter what.
1367 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1370 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1373 /* OK. We couldn't find a significantly lower energy.
1374 * If beta==0 this was steepest descent, and then we give up.
1375 * If not, set beta=0 and restart with steepest descent before quitting.
1385 /* Reset memory before giving up */
1391 /* Select min energy state of A & C, put the best in B.
1393 if (s_c->epot < s_a->epot)
1397 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1398 s_c->epot, s_a->epot);
1400 swap_em_state(s_b, s_c);
1408 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1409 s_a->epot, s_c->epot);
1411 swap_em_state(s_b, s_a);
1421 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1424 swap_em_state(s_b, s_c);
1429 /* new search direction */
1430 /* beta = 0 means forget all memory and restart with steepest descents. */
1431 if (nstcg && ((step % nstcg) == 0))
1437 /* s_min->fnorm cannot be zero, because then we would have converged
1441 /* Polak-Ribiere update.
1442 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1444 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1446 /* Limit beta to prevent oscillations */
1447 if (fabs(beta) > 5.0)
1453 /* update positions */
1454 swap_em_state(s_min, s_b);
1457 /* Print it if necessary */
1462 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1463 step, s_min->epot, s_min->fnorm/sqrt(state_global->natoms),
1464 s_min->fmax, s_min->a_fmax+1);
1466 /* Store the new (lower) energies */
1467 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1468 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1469 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1471 do_log = do_per_step(step, inputrec->nstlog);
1472 do_ene = do_per_step(step, inputrec->nstenergy);
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 /* Stop when the maximum force lies below tolerance.
1483 * If we have reached machine precision, converged is already set to true.
1485 converged = converged || (s_min->fmax < inputrec->em_tol);
1487 } /* End of the loop */
1491 step--; /* we never took that last step in this case */
1494 if (s_min->fmax > inputrec->em_tol)
1498 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1499 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1506 /* If we printed energy and/or logfile last step (which was the last step)
1507 * we don't have to do it again, but otherwise print the final values.
1511 /* Write final value to log since we didn't do anything the last step */
1512 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1514 if (!do_ene || !do_log)
1516 /* Write final energy file entries */
1517 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1518 !do_log ? fplog : NULL, step, step, eprNORMAL,
1519 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1523 /* Print some stuff... */
1526 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1530 * For accurate normal mode calculation it is imperative that we
1531 * store the last conformation into the full precision binary trajectory.
1533 * However, we should only do it if we did NOT already write this step
1534 * above (which we did if do_x or do_f was true).
1536 do_x = !do_per_step(step, inputrec->nstxout);
1537 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1539 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1540 top_global, inputrec, step,
1541 s_min, state_global, f_global);
1543 fnormn = s_min->fnorm/sqrt(state_global->natoms);
1547 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1548 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1549 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1550 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1552 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1555 finish_em(cr, outf, walltime_accounting, wcycle);
1557 /* To print the actual number of steps we needed somewhere */
1558 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1561 } /* That's all folks */
1564 double do_lbfgs(FILE *fplog, t_commrec *cr,
1565 int nfile, const t_filenm fnm[],
1566 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
1567 int gmx_unused nstglobalcomm,
1568 gmx_vsite_t *vsite, gmx_constr_t constr,
1569 int gmx_unused stepout,
1570 t_inputrec *inputrec,
1571 gmx_mtop_t *top_global, t_fcdata *fcd,
1574 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1575 gmx_edsam_t gmx_unused ed,
1577 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1578 gmx_membed_t gmx_unused membed,
1579 real gmx_unused cpt_period, real gmx_unused max_hours,
1580 const char gmx_unused *deviceOptions,
1581 unsigned long gmx_unused Flags,
1582 gmx_walltime_accounting_t walltime_accounting)
1584 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1586 gmx_localtop_t *top;
1587 gmx_enerdata_t *enerd;
1589 gmx_global_stat_t gstat;
1592 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1593 double stepsize, gpa, gpb, gpc, tmp, minstep;
1594 real *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1595 real *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1596 real a, b, c, maxdelta, delta;
1597 real diag, Epot0, Epot, EpotA, EpotB, EpotC;
1598 real dgdx, dgdg, sq, yr, beta;
1600 gmx_bool converged, first;
1603 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1605 int start, end, number_steps;
1607 int i, k, m, n, nfmax, gf, step;
1614 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1619 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).");
1622 n = 3*state->natoms;
1623 nmaxcorr = inputrec->nbfgscorr;
1625 /* Allocate memory */
1626 /* Use pointers to real so we dont have to loop over both atoms and
1627 * dimensions all the time...
1628 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1629 * that point to the same memory.
1642 snew(rho, nmaxcorr);
1643 snew(alpha, nmaxcorr);
1646 for (i = 0; i < nmaxcorr; i++)
1652 for (i = 0; i < nmaxcorr; i++)
1661 init_em(fplog, LBFGS, cr, inputrec,
1662 state, top_global, &ems, &top, &f, &f_global,
1663 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1664 nfile, fnm, &outf, &mdebin);
1665 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1666 * so we free some memory again.
1671 xx = (real *)state->x;
1674 start = mdatoms->start;
1675 end = mdatoms->homenr + start;
1677 /* Print to log file */
1678 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1680 do_log = do_ene = do_x = do_f = TRUE;
1682 /* Max number of steps */
1683 number_steps = inputrec->nsteps;
1685 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1687 for (i = start; i < end; i++)
1689 if (mdatoms->cFREEZE)
1691 gf = mdatoms->cFREEZE[i];
1693 for (m = 0; m < DIM; m++)
1695 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1700 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1704 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1709 construct_vsites(vsite, state->x, 1, NULL,
1710 top->idef.iparams, top->idef.il,
1711 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1714 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1715 /* do_force always puts the charge groups in the box and shifts again
1716 * We do not unshift, so molecules are always whole
1721 evaluate_energy(fplog, cr,
1722 top_global, &ems, top,
1723 inputrec, nrnb, wcycle, gstat,
1724 vsite, constr, fcd, graph, mdatoms, fr,
1725 mu_tot, enerd, vir, pres, -1, TRUE);
1730 /* Copy stuff to the energy bin for easy printing etc. */
1731 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1732 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
1733 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1735 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
1736 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1737 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1741 /* This is the starting energy */
1742 Epot = enerd->term[F_EPOT];
1748 /* Set the initial step.
1749 * since it will be multiplied by the non-normalized search direction
1750 * vector (force vector the first time), we scale it by the
1751 * norm of the force.
1756 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1757 fprintf(stderr, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1758 fprintf(stderr, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1759 fprintf(stderr, "\n");
1760 /* and copy to the log file too... */
1761 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1762 fprintf(fplog, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1763 fprintf(fplog, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1764 fprintf(fplog, "\n");
1768 for (i = 0; i < n; i++)
1772 dx[point][i] = ff[i]; /* Initial search direction */
1780 stepsize = 1.0/fnorm;
1783 /* Start the loop over BFGS steps.
1784 * Each successful step is counted, and we continue until
1785 * we either converge or reach the max number of steps.
1790 /* Set the gradient from the force */
1792 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1795 /* Write coordinates if necessary */
1796 do_x = do_per_step(step, inputrec->nstxout);
1797 do_f = do_per_step(step, inputrec->nstfout);
1802 mdof_flags |= MDOF_X;
1807 mdof_flags |= MDOF_F;
1810 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1811 top_global, step, (real)step, state, state, f, f);
1813 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1815 /* pointer to current direction - point=0 first time here */
1818 /* calculate line gradient */
1819 for (gpa = 0, i = 0; i < n; i++)
1824 /* Calculate minimum allowed stepsize, before the average (norm)
1825 * relative change in coordinate is smaller than precision
1827 for (minstep = 0, i = 0; i < n; i++)
1837 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1839 if (stepsize < minstep)
1845 /* Store old forces and coordinates */
1846 for (i = 0; i < n; i++)
1855 for (i = 0; i < n; i++)
1860 /* Take a step downhill.
1861 * In theory, we should minimize the function along this direction.
1862 * That is quite possible, but it turns out to take 5-10 function evaluations
1863 * for each line. However, we dont really need to find the exact minimum -
1864 * it is much better to start a new BFGS step in a modified direction as soon
1865 * as we are close to it. This will save a lot of energy evaluations.
1867 * In practice, we just try to take a single step.
1868 * If it worked (i.e. lowered the energy), we increase the stepsize but
1869 * the continue straight to the next BFGS step without trying to find any minimum.
1870 * If it didn't work (higher energy), there must be a minimum somewhere between
1871 * the old position and the new one.
1873 * Due to the finite numerical accuracy, it turns out that it is a good idea
1874 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1875 * This leads to lower final energies in the tests I've done. / Erik
1880 c = a + stepsize; /* reference position along line is zero */
1882 /* Check stepsize first. We do not allow displacements
1883 * larger than emstep.
1889 for (i = 0; i < n; i++)
1892 if (delta > maxdelta)
1897 if (maxdelta > inputrec->em_stepsize)
1902 while (maxdelta > inputrec->em_stepsize);
1904 /* Take a trial step */
1905 for (i = 0; i < n; i++)
1907 xc[i] = lastx[i] + c*s[i];
1911 /* Calculate energy for the trial step */
1912 ems.s.x = (rvec *)xc;
1914 evaluate_energy(fplog, cr,
1915 top_global, &ems, top,
1916 inputrec, nrnb, wcycle, gstat,
1917 vsite, constr, fcd, graph, mdatoms, fr,
1918 mu_tot, enerd, vir, pres, step, FALSE);
1921 /* Calc derivative along line */
1922 for (gpc = 0, i = 0; i < n; i++)
1924 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1926 /* Sum the gradient along the line across CPUs */
1929 gmx_sumd(1, &gpc, cr);
1932 /* This is the max amount of increase in energy we tolerate */
1933 tmp = sqrt(GMX_REAL_EPS)*fabs(EpotA);
1935 /* Accept the step if the energy is lower, or if it is not significantly higher
1936 * and the line derivative is still negative.
1938 if (EpotC < EpotA || (gpc < 0 && EpotC < (EpotA+tmp)))
1941 /* Great, we found a better energy. Increase step for next iteration
1942 * if we are still going down, decrease it otherwise
1946 stepsize *= 1.618034; /* The golden section */
1950 stepsize *= 0.618034; /* 1/golden section */
1955 /* New energy is the same or higher. We will have to do some work
1956 * to find a smaller value in the interval. Take smaller step next time!
1959 stepsize *= 0.618034;
1962 /* OK, if we didn't find a lower value we will have to locate one now - there must
1963 * be one in the interval [a=0,c].
1964 * The same thing is valid here, though: Don't spend dozens of iterations to find
1965 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1966 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1968 * I also have a safeguard for potentially really patological functions so we never
1969 * take more than 20 steps before we give up ...
1971 * If we already found a lower value we just skip this step and continue to the update.
1980 /* Select a new trial point.
1981 * If the derivatives at points a & c have different sign we interpolate to zero,
1982 * otherwise just do a bisection.
1985 if (gpa < 0 && gpc > 0)
1987 b = a + gpa*(a-c)/(gpc-gpa);
1994 /* safeguard if interpolation close to machine accuracy causes errors:
1995 * never go outside the interval
1997 if (b <= a || b >= c)
2002 /* Take a trial step */
2003 for (i = 0; i < n; i++)
2005 xb[i] = lastx[i] + b*s[i];
2009 /* Calculate energy for the trial step */
2010 ems.s.x = (rvec *)xb;
2012 evaluate_energy(fplog, cr,
2013 top_global, &ems, top,
2014 inputrec, nrnb, wcycle, gstat,
2015 vsite, constr, fcd, graph, mdatoms, fr,
2016 mu_tot, enerd, vir, pres, step, FALSE);
2021 for (gpb = 0, i = 0; i < n; i++)
2023 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2026 /* Sum the gradient along the line across CPUs */
2029 gmx_sumd(1, &gpb, cr);
2032 /* Keep one of the intervals based on the value of the derivative at the new point */
2035 /* Replace c endpoint with b */
2039 /* swap coord pointers b/c */
2049 /* Replace a endpoint with b */
2053 /* swap coord pointers a/b */
2063 * Stop search as soon as we find a value smaller than the endpoints,
2064 * or if the tolerance is below machine precision.
2065 * Never run more than 20 steps, no matter what.
2069 while ((EpotB > EpotA || EpotB > EpotC) && (nminstep < 20));
2071 if (fabs(EpotB-Epot0) < GMX_REAL_EPS || nminstep >= 20)
2073 /* OK. We couldn't find a significantly lower energy.
2074 * If ncorr==0 this was steepest descent, and then we give up.
2075 * If not, reset memory to restart as steepest descent before quitting.
2087 /* Search in gradient direction */
2088 for (i = 0; i < n; i++)
2090 dx[point][i] = ff[i];
2092 /* Reset stepsize */
2093 stepsize = 1.0/fnorm;
2098 /* Select min energy state of A & C, put the best in xx/ff/Epot
2104 for (i = 0; i < n; i++)
2115 for (i = 0; i < n; i++)
2129 for (i = 0; i < n; i++)
2137 /* Update the memory information, and calculate a new
2138 * approximation of the inverse hessian
2141 /* Have new data in Epot, xx, ff */
2142 if (ncorr < nmaxcorr)
2147 for (i = 0; i < n; i++)
2149 dg[point][i] = lastf[i]-ff[i];
2150 dx[point][i] *= stepsize;
2155 for (i = 0; i < n; i++)
2157 dgdg += dg[point][i]*dg[point][i];
2158 dgdx += dg[point][i]*dx[point][i];
2163 rho[point] = 1.0/dgdx;
2166 if (point >= nmaxcorr)
2172 for (i = 0; i < n; i++)
2179 /* Recursive update. First go back over the memory points */
2180 for (k = 0; k < ncorr; k++)
2189 for (i = 0; i < n; i++)
2191 sq += dx[cp][i]*p[i];
2194 alpha[cp] = rho[cp]*sq;
2196 for (i = 0; i < n; i++)
2198 p[i] -= alpha[cp]*dg[cp][i];
2202 for (i = 0; i < n; i++)
2207 /* And then go forward again */
2208 for (k = 0; k < ncorr; k++)
2211 for (i = 0; i < n; i++)
2213 yr += p[i]*dg[cp][i];
2217 beta = alpha[cp]-beta;
2219 for (i = 0; i < n; i++)
2221 p[i] += beta*dx[cp][i];
2231 for (i = 0; i < n; i++)
2235 dx[point][i] = p[i];
2245 /* Test whether the convergence criterion is met */
2246 get_f_norm_max(cr, &(inputrec->opts), mdatoms, f, &fnorm, &fmax, &nfmax);
2248 /* Print it if necessary */
2253 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2254 step, Epot, fnorm/sqrt(state->natoms), fmax, nfmax+1);
2256 /* Store the new (lower) energies */
2257 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2258 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
2259 NULL, NULL, vir, pres, NULL, mu_tot, constr);
2260 do_log = do_per_step(step, inputrec->nstlog);
2261 do_ene = do_per_step(step, inputrec->nstenergy);
2264 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2266 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2267 do_log ? fplog : NULL, step, step, eprNORMAL,
2268 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2271 /* Stop when the maximum force lies below tolerance.
2272 * If we have reached machine precision, converged is already set to true.
2275 converged = converged || (fmax < inputrec->em_tol);
2277 } /* End of the loop */
2281 step--; /* we never took that last step in this case */
2284 if (fmax > inputrec->em_tol)
2288 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2289 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2294 /* If we printed energy and/or logfile last step (which was the last step)
2295 * we don't have to do it again, but otherwise print the final values.
2297 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2299 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2301 if (!do_ene || !do_log) /* Write final energy file entries */
2303 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2304 !do_log ? fplog : NULL, step, step, eprNORMAL,
2305 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2308 /* Print some stuff... */
2311 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2315 * For accurate normal mode calculation it is imperative that we
2316 * store the last conformation into the full precision binary trajectory.
2318 * However, we should only do it if we did NOT already write this step
2319 * above (which we did if do_x or do_f was true).
2321 do_x = !do_per_step(step, inputrec->nstxout);
2322 do_f = !do_per_step(step, inputrec->nstfout);
2323 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2324 top_global, inputrec, step,
2329 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2330 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2331 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2332 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2334 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2337 finish_em(cr, outf, walltime_accounting, wcycle);
2339 /* To print the actual number of steps we needed somewhere */
2340 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2343 } /* That's all folks */
2346 double do_steep(FILE *fplog, t_commrec *cr,
2347 int nfile, const t_filenm fnm[],
2348 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2349 int gmx_unused nstglobalcomm,
2350 gmx_vsite_t *vsite, gmx_constr_t constr,
2351 int gmx_unused stepout,
2352 t_inputrec *inputrec,
2353 gmx_mtop_t *top_global, t_fcdata *fcd,
2354 t_state *state_global,
2356 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2357 gmx_edsam_t gmx_unused ed,
2359 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2360 gmx_membed_t gmx_unused membed,
2361 real gmx_unused cpt_period, real gmx_unused max_hours,
2362 const char gmx_unused *deviceOptions,
2363 unsigned long gmx_unused Flags,
2364 gmx_walltime_accounting_t walltime_accounting)
2366 const char *SD = "Steepest Descents";
2367 em_state_t *s_min, *s_try;
2369 gmx_localtop_t *top;
2370 gmx_enerdata_t *enerd;
2372 gmx_global_stat_t gstat;
2374 real stepsize, constepsize;
2378 gmx_bool bDone, bAbort, do_x, do_f;
2383 int steps_accepted = 0;
2387 s_min = init_em_state();
2388 s_try = init_em_state();
2390 /* Init em and store the local state in s_try */
2391 init_em(fplog, SD, cr, inputrec,
2392 state_global, top_global, s_try, &top, &f, &f_global,
2393 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2394 nfile, fnm, &outf, &mdebin);
2396 /* Print to log file */
2397 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2399 /* Set variables for stepsize (in nm). This is the largest
2400 * step that we are going to make in any direction.
2402 ustep = inputrec->em_stepsize;
2405 /* Max number of steps */
2406 nsteps = inputrec->nsteps;
2410 /* Print to the screen */
2411 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2415 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2418 /**** HERE STARTS THE LOOP ****
2419 * count is the counter for the number of steps
2420 * bDone will be TRUE when the minimization has converged
2421 * bAbort will be TRUE when nsteps steps have been performed or when
2422 * the stepsize becomes smaller than is reasonable for machine precision
2427 while (!bDone && !bAbort)
2429 bAbort = (nsteps >= 0) && (count == nsteps);
2431 /* set new coordinates, except for first step */
2434 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2435 s_min, stepsize, s_min->f, s_try,
2436 constr, top, nrnb, wcycle, count);
2439 evaluate_energy(fplog, cr,
2440 top_global, s_try, top,
2441 inputrec, nrnb, wcycle, gstat,
2442 vsite, constr, fcd, graph, mdatoms, fr,
2443 mu_tot, enerd, vir, pres, count, count == 0);
2447 print_ebin_header(fplog, count, count, s_try->s.lambda[efptFEP]);
2452 s_min->epot = s_try->epot + 1;
2455 /* Print it if necessary */
2460 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2461 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2462 (s_try->epot < s_min->epot) ? '\n' : '\r');
2465 if (s_try->epot < s_min->epot)
2467 /* Store the new (lower) energies */
2468 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2469 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2470 s_try->s.box, NULL, NULL, vir, pres, NULL, mu_tot, constr);
2471 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2472 do_per_step(steps_accepted, inputrec->nstdisreout),
2473 do_per_step(steps_accepted, inputrec->nstorireout),
2474 fplog, count, count, eprNORMAL, TRUE,
2475 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2480 /* Now if the new energy is smaller than the previous...
2481 * or if this is the first step!
2482 * or if we did random steps!
2485 if ( (count == 0) || (s_try->epot < s_min->epot) )
2489 /* Test whether the convergence criterion is met... */
2490 bDone = (s_try->fmax < inputrec->em_tol);
2492 /* Copy the arrays for force, positions and energy */
2493 /* The 'Min' array always holds the coords and forces of the minimal
2495 swap_em_state(s_min, s_try);
2501 /* Write to trn, if necessary */
2502 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2503 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2504 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
2505 top_global, inputrec, count,
2506 s_min, state_global, f_global);
2510 /* If energy is not smaller make the step smaller... */
2513 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2515 /* Reload the old state */
2516 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2517 s_min, top, mdatoms, fr, vsite, constr,
2522 /* Determine new step */
2523 stepsize = ustep/s_min->fmax;
2525 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2527 if (count == nsteps || ustep < 1e-12)
2529 if (count == nsteps || ustep < 1e-6)
2534 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2535 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2541 } /* End of the loop */
2543 /* Print some shit... */
2546 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2548 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2549 top_global, inputrec, count,
2550 s_min, state_global, f_global);
2552 fnormn = s_min->fnorm/sqrt(state_global->natoms);
2556 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2557 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2558 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2559 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2562 finish_em(cr, outf, walltime_accounting, wcycle);
2564 /* To print the actual number of steps we needed somewhere */
2565 inputrec->nsteps = count;
2567 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2570 } /* That's all folks */
2573 double do_nm(FILE *fplog, t_commrec *cr,
2574 int nfile, const t_filenm fnm[],
2575 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2576 int gmx_unused nstglobalcomm,
2577 gmx_vsite_t *vsite, gmx_constr_t constr,
2578 int gmx_unused stepout,
2579 t_inputrec *inputrec,
2580 gmx_mtop_t *top_global, t_fcdata *fcd,
2581 t_state *state_global,
2583 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2584 gmx_edsam_t gmx_unused ed,
2586 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2587 gmx_membed_t gmx_unused membed,
2588 real gmx_unused cpt_period, real gmx_unused max_hours,
2589 const char gmx_unused *deviceOptions,
2590 unsigned long gmx_unused Flags,
2591 gmx_walltime_accounting_t walltime_accounting)
2593 const char *NM = "Normal Mode Analysis";
2595 int natoms, atom, d;
2598 gmx_localtop_t *top;
2599 gmx_enerdata_t *enerd;
2601 gmx_global_stat_t gstat;
2603 real t, t0, lambda, lam0;
2608 gmx_bool bSparse; /* use sparse matrix storage format */
2610 gmx_sparsematrix_t * sparse_matrix = NULL;
2611 real * full_matrix = NULL;
2612 em_state_t * state_work;
2614 /* added with respect to mdrun */
2615 int i, j, k, row, col;
2616 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2622 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2625 state_work = init_em_state();
2627 /* Init em and store the local state in state_minimum */
2628 init_em(fplog, NM, cr, inputrec,
2629 state_global, top_global, state_work, &top,
2631 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2632 nfile, fnm, &outf, NULL);
2634 natoms = top_global->natoms;
2642 "NOTE: This version of Gromacs has been compiled in single precision,\n"
2643 " which MIGHT not be accurate enough for normal mode analysis.\n"
2644 " Gromacs now uses sparse matrix storage, so the memory requirements\n"
2645 " are fairly modest even if you recompile in double precision.\n\n");
2649 /* Check if we can/should use sparse storage format.
2651 * Sparse format is only useful when the Hessian itself is sparse, which it
2652 * will be when we use a cutoff.
2653 * For small systems (n<1000) it is easier to always use full matrix format, though.
2655 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2657 md_print_info(cr, fplog, "Non-cutoff electrostatics used, forcing full Hessian format.\n");
2660 else if (top_global->natoms < 1000)
2662 md_print_info(cr, fplog, "Small system size (N=%d), using full Hessian format.\n", top_global->natoms);
2667 md_print_info(cr, fplog, "Using compressed symmetric sparse Hessian format.\n");
2673 sz = DIM*top_global->natoms;
2675 fprintf(stderr, "Allocating Hessian memory...\n\n");
2679 sparse_matrix = gmx_sparsematrix_init(sz);
2680 sparse_matrix->compressed_symmetric = TRUE;
2684 snew(full_matrix, sz*sz);
2688 /* Initial values */
2689 t0 = inputrec->init_t;
2690 lam0 = inputrec->fepvals->init_lambda;
2698 /* Write start time and temperature */
2699 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2701 /* fudge nr of steps to nr of atoms */
2702 inputrec->nsteps = natoms*2;
2706 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2707 *(top_global->name), (int)inputrec->nsteps);
2710 nnodes = cr->nnodes;
2712 /* Make evaluate_energy do a single node force calculation */
2714 evaluate_energy(fplog, cr,
2715 top_global, state_work, top,
2716 inputrec, nrnb, wcycle, gstat,
2717 vsite, constr, fcd, graph, mdatoms, fr,
2718 mu_tot, enerd, vir, pres, -1, TRUE);
2719 cr->nnodes = nnodes;
2721 /* if forces are not small, warn user */
2722 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, state_work);
2724 md_print_info(cr, fplog, "Maximum force:%12.5e\n", state_work->fmax);
2725 if (state_work->fmax > 1.0e-3)
2727 md_print_info(cr, fplog,
2728 "The force is probably not small enough to "
2729 "ensure that you are at a minimum.\n"
2730 "Be aware that negative eigenvalues may occur\n"
2731 "when the resulting matrix is diagonalized.\n\n");
2734 /***********************************************************
2736 * Loop over all pairs in matrix
2738 * do_force called twice. Once with positive and
2739 * once with negative displacement
2741 ************************************************************/
2743 /* Steps are divided one by one over the nodes */
2744 for (atom = cr->nodeid; atom < natoms; atom += nnodes)
2747 for (d = 0; d < DIM; d++)
2749 x_min = state_work->s.x[atom][d];
2751 state_work->s.x[atom][d] = x_min - der_range;
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, atom*2, FALSE);
2761 for (i = 0; i < natoms; i++)
2763 copy_rvec(state_work->f[i], fneg[i]);
2766 state_work->s.x[atom][d] = x_min + der_range;
2768 evaluate_energy(fplog, cr,
2769 top_global, state_work, top,
2770 inputrec, nrnb, wcycle, gstat,
2771 vsite, constr, fcd, graph, mdatoms, fr,
2772 mu_tot, enerd, vir, pres, atom*2+1, FALSE);
2773 cr->nnodes = nnodes;
2775 /* x is restored to original */
2776 state_work->s.x[atom][d] = x_min;
2778 for (j = 0; j < natoms; j++)
2780 for (k = 0; (k < DIM); k++)
2783 -(state_work->f[j][k] - fneg[j][k])/(2*der_range);
2791 #define mpi_type MPI_DOUBLE
2793 #define mpi_type MPI_FLOAT
2795 MPI_Send(dfdx[0], natoms*DIM, mpi_type, MASTERNODE(cr), cr->nodeid,
2796 cr->mpi_comm_mygroup);
2801 for (node = 0; (node < nnodes && atom+node < natoms); node++)
2807 MPI_Recv(dfdx[0], natoms*DIM, mpi_type, node, node,
2808 cr->mpi_comm_mygroup, &stat);
2813 row = (atom + node)*DIM + d;
2815 for (j = 0; j < natoms; j++)
2817 for (k = 0; k < DIM; k++)
2823 if (col >= row && dfdx[j][k] != 0.0)
2825 gmx_sparsematrix_increment_value(sparse_matrix,
2826 row, col, dfdx[j][k]);
2831 full_matrix[row*sz+col] = dfdx[j][k];
2838 if (bVerbose && fplog)
2843 /* write progress */
2844 if (MASTER(cr) && bVerbose)
2846 fprintf(stderr, "\rFinished step %d out of %d",
2847 min(atom+nnodes, natoms), natoms);
2854 fprintf(stderr, "\n\nWriting Hessian...\n");
2855 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2858 finish_em(cr, outf, walltime_accounting, wcycle);
2860 walltime_accounting_set_nsteps_done(walltime_accounting, natoms*2);