When writing TNG include file closing in wallcycle.
[gromacs.git] / src / gromacs / mdlib / minimize.c
blob3cc1bbc8a08c6a9d65f6d324369dc49cc6addcd3
1 /*
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.
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
41 #include <string.h>
42 #include <time.h>
43 #include <math.h>
44 #include "sysstuff.h"
45 #include "gromacs/utility/cstringutil.h"
46 #include "network.h"
47 #include "gromacs/utility/smalloc.h"
48 #include "nrnb.h"
49 #include "main.h"
50 #include "force.h"
51 #include "macros.h"
52 #include "names.h"
53 #include "gmx_fatal.h"
54 #include "txtdump.h"
55 #include "typedefs.h"
56 #include "update.h"
57 #include "constr.h"
58 #include "vec.h"
59 #include "tgroup.h"
60 #include "mdebin.h"
61 #include "vsite.h"
62 #include "force.h"
63 #include "mdrun.h"
64 #include "md_support.h"
65 #include "sim_util.h"
66 #include "domdec.h"
67 #include "mdatoms.h"
68 #include "ns.h"
69 #include "mtop_util.h"
70 #include "pme.h"
71 #include "bondf.h"
72 #include "gmx_omp_nthreads.h"
73 #include "md_logging.h"
75 #include "gromacs/fileio/confio.h"
76 #include "gromacs/fileio/trajectory_writing.h"
77 #include "gromacs/linearalgebra/mtxio.h"
78 #include "gromacs/linearalgebra/sparsematrix.h"
79 #include "gromacs/timing/wallcycle.h"
80 #include "gromacs/timing/walltime_accounting.h"
81 #include "gromacs/imd/imd.h"
83 typedef struct {
84 t_state s;
85 rvec *f;
86 real epot;
87 real fnorm;
88 real fmax;
89 int a_fmax;
90 } em_state_t;
92 static em_state_t *init_em_state()
94 em_state_t *ems;
96 snew(ems, 1);
98 /* does this need to be here? Should the array be declared differently (staticaly)in the state definition? */
99 snew(ems->s.lambda, efptNR);
101 return ems;
104 static void print_em_start(FILE *fplog,
105 t_commrec *cr,
106 gmx_walltime_accounting_t walltime_accounting,
107 gmx_wallcycle_t wcycle,
108 const char *name)
110 walltime_accounting_start(walltime_accounting);
111 wallcycle_start(wcycle, ewcRUN);
112 print_start(fplog, cr, walltime_accounting, name);
114 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
115 gmx_wallcycle_t wcycle)
117 wallcycle_stop(wcycle, ewcRUN);
119 walltime_accounting_end(walltime_accounting);
122 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
124 fprintf(out, "\n");
125 fprintf(out, "%s:\n", minimizer);
126 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
127 fprintf(out, " Number of steps = %12d\n", nsteps);
130 static void warn_step(FILE *fp, real ftol, gmx_bool bLastStep, gmx_bool bConstrain)
132 char buffer[2048];
133 if (bLastStep)
135 sprintf(buffer,
136 "\nEnergy minimization reached the maximum number "
137 "of steps before the forces reached the requested "
138 "precision Fmax < %g.\n", ftol);
140 else
142 sprintf(buffer,
143 "\nEnergy minimization has stopped, but the forces have "
144 "not converged to the requested precision Fmax < %g (which "
145 "may not be possible for your system). It stopped "
146 "because the algorithm tried to make a new step whose size "
147 "was too small, or there was no change in the energy since "
148 "last step. Either way, we regard the minimization as "
149 "converged to within the available machine precision, "
150 "given your starting configuration and EM parameters.\n%s%s",
151 ftol,
152 sizeof(real) < sizeof(double) ?
153 "\nDouble precision normally gives you higher accuracy, but "
154 "this is often not needed for preparing to run molecular "
155 "dynamics.\n" :
157 bConstrain ?
158 "You might need to increase your constraint accuracy, or turn\n"
159 "off constraints altogether (set constraints = none in mdp file)\n" :
160 "");
162 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
167 static void print_converged(FILE *fp, const char *alg, real ftol,
168 gmx_int64_t count, gmx_bool bDone, gmx_int64_t nsteps,
169 real epot, real fmax, int nfmax, real fnorm)
171 char buf[STEPSTRSIZE];
173 if (bDone)
175 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
176 alg, ftol, gmx_step_str(count, buf));
178 else if (count < nsteps)
180 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
181 "but did not reach the requested Fmax < %g.\n",
182 alg, gmx_step_str(count, buf), ftol);
184 else
186 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
187 alg, ftol, gmx_step_str(count, buf));
190 #ifdef GMX_DOUBLE
191 fprintf(fp, "Potential Energy = %21.14e\n", epot);
192 fprintf(fp, "Maximum force = %21.14e on atom %d\n", fmax, nfmax+1);
193 fprintf(fp, "Norm of force = %21.14e\n", fnorm);
194 #else
195 fprintf(fp, "Potential Energy = %14.7e\n", epot);
196 fprintf(fp, "Maximum force = %14.7e on atom %d\n", fmax, nfmax+1);
197 fprintf(fp, "Norm of force = %14.7e\n", fnorm);
198 #endif
201 static void get_f_norm_max(t_commrec *cr,
202 t_grpopts *opts, t_mdatoms *mdatoms, rvec *f,
203 real *fnorm, real *fmax, int *a_fmax)
205 double fnorm2, *sum;
206 real fmax2, fmax2_0, fam;
207 int la_max, a_max, start, end, i, m, gf;
209 /* This routine finds the largest force and returns it.
210 * On parallel machines the global max is taken.
212 fnorm2 = 0;
213 fmax2 = 0;
214 la_max = -1;
215 gf = 0;
216 start = 0;
217 end = mdatoms->homenr;
218 if (mdatoms->cFREEZE)
220 for (i = start; i < end; i++)
222 gf = mdatoms->cFREEZE[i];
223 fam = 0;
224 for (m = 0; m < DIM; m++)
226 if (!opts->nFreeze[gf][m])
228 fam += sqr(f[i][m]);
231 fnorm2 += fam;
232 if (fam > fmax2)
234 fmax2 = fam;
235 la_max = i;
239 else
241 for (i = start; i < end; i++)
243 fam = norm2(f[i]);
244 fnorm2 += fam;
245 if (fam > fmax2)
247 fmax2 = fam;
248 la_max = i;
253 if (la_max >= 0 && DOMAINDECOMP(cr))
255 a_max = cr->dd->gatindex[la_max];
257 else
259 a_max = la_max;
261 if (PAR(cr))
263 snew(sum, 2*cr->nnodes+1);
264 sum[2*cr->nodeid] = fmax2;
265 sum[2*cr->nodeid+1] = a_max;
266 sum[2*cr->nnodes] = fnorm2;
267 gmx_sumd(2*cr->nnodes+1, sum, cr);
268 fnorm2 = sum[2*cr->nnodes];
269 /* Determine the global maximum */
270 for (i = 0; i < cr->nnodes; i++)
272 if (sum[2*i] > fmax2)
274 fmax2 = sum[2*i];
275 a_max = (int)(sum[2*i+1] + 0.5);
278 sfree(sum);
281 if (fnorm)
283 *fnorm = sqrt(fnorm2);
285 if (fmax)
287 *fmax = sqrt(fmax2);
289 if (a_fmax)
291 *a_fmax = a_max;
295 static void get_state_f_norm_max(t_commrec *cr,
296 t_grpopts *opts, t_mdatoms *mdatoms,
297 em_state_t *ems)
299 get_f_norm_max(cr, opts, mdatoms, ems->f, &ems->fnorm, &ems->fmax, &ems->a_fmax);
302 void init_em(FILE *fplog, const char *title,
303 t_commrec *cr, t_inputrec *ir,
304 t_state *state_global, gmx_mtop_t *top_global,
305 em_state_t *ems, gmx_localtop_t **top,
306 rvec **f, rvec **f_global,
307 t_nrnb *nrnb, rvec mu_tot,
308 t_forcerec *fr, gmx_enerdata_t **enerd,
309 t_graph **graph, t_mdatoms *mdatoms, gmx_global_stat_t *gstat,
310 gmx_vsite_t *vsite, gmx_constr_t constr,
311 int nfile, const t_filenm fnm[],
312 gmx_mdoutf_t *outf, t_mdebin **mdebin,
313 int imdport, unsigned long gmx_unused Flags,
314 gmx_wallcycle_t wcycle)
316 int i;
317 real dvdl_constr;
319 if (fplog)
321 fprintf(fplog, "Initiating %s\n", title);
324 state_global->ngtc = 0;
326 /* Initialize lambda variables */
327 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, NULL);
329 init_nrnb(nrnb);
331 /* Interactive molecular dynamics */
332 init_IMD(ir, cr, top_global, fplog, 1, state_global->x,
333 nfile, fnm, NULL, imdport, Flags);
335 if (DOMAINDECOMP(cr))
337 *top = dd_init_local_top(top_global);
339 dd_init_local_state(cr->dd, state_global, &ems->s);
341 *f = NULL;
343 /* Distribute the charge groups over the nodes from the master node */
344 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
345 state_global, top_global, ir,
346 &ems->s, &ems->f, mdatoms, *top,
347 fr, vsite, NULL, constr,
348 nrnb, NULL, FALSE);
349 dd_store_state(cr->dd, &ems->s);
351 if (ir->nstfout)
353 snew(*f_global, top_global->natoms);
355 else
357 *f_global = NULL;
359 *graph = NULL;
361 else
363 snew(*f, top_global->natoms);
365 /* Just copy the state */
366 ems->s = *state_global;
367 snew(ems->s.x, ems->s.nalloc);
368 snew(ems->f, ems->s.nalloc);
369 for (i = 0; i < state_global->natoms; i++)
371 copy_rvec(state_global->x[i], ems->s.x[i]);
373 copy_mat(state_global->box, ems->s.box);
375 *top = gmx_mtop_generate_local_top(top_global, ir);
376 *f_global = *f;
378 forcerec_set_excl_load(fr, *top);
380 setup_bonded_threading(fr, &(*top)->idef);
382 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
384 *graph = mk_graph(fplog, &((*top)->idef), 0, top_global->natoms, FALSE, FALSE);
386 else
388 *graph = NULL;
391 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
392 update_mdatoms(mdatoms, state_global->lambda[efptFEP]);
394 if (vsite)
396 set_vsite_top(vsite, *top, mdatoms, cr);
400 if (constr)
402 if (ir->eConstrAlg == econtSHAKE &&
403 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
405 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
406 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
409 if (!DOMAINDECOMP(cr))
411 set_constraints(constr, *top, ir, mdatoms, cr);
414 if (!ir->bContinuation)
416 /* Constrain the starting coordinates */
417 dvdl_constr = 0;
418 constrain(PAR(cr) ? NULL : fplog, TRUE, TRUE, constr, &(*top)->idef,
419 ir, NULL, cr, -1, 0, 1.0, mdatoms,
420 ems->s.x, ems->s.x, NULL, fr->bMolPBC, ems->s.box,
421 ems->s.lambda[efptFEP], &dvdl_constr,
422 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
426 if (PAR(cr))
428 *gstat = global_stat_init(ir);
431 *outf = init_mdoutf(fplog, nfile, fnm, 0, cr, ir, top_global, NULL, wcycle);
433 snew(*enerd, 1);
434 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
435 *enerd);
437 if (mdebin != NULL)
439 /* Init bin for energy stuff */
440 *mdebin = init_mdebin(mdoutf_get_fp_ene(*outf), top_global, ir, NULL);
443 clear_rvec(mu_tot);
444 calc_shifts(ems->s.box, fr->shift_vec);
447 static void finish_em(t_commrec *cr, gmx_mdoutf_t outf,
448 gmx_walltime_accounting_t walltime_accounting,
449 gmx_wallcycle_t wcycle)
451 if (!(cr->duty & DUTY_PME))
453 /* Tell the PME only node to finish */
454 gmx_pme_send_finish(cr);
457 done_mdoutf(outf);
459 em_time_end(walltime_accounting, wcycle);
462 static void swap_em_state(em_state_t *ems1, em_state_t *ems2)
464 em_state_t tmp;
466 tmp = *ems1;
467 *ems1 = *ems2;
468 *ems2 = tmp;
471 static void copy_em_coords(em_state_t *ems, t_state *state)
473 int i;
475 for (i = 0; (i < state->natoms); i++)
477 copy_rvec(ems->s.x[i], state->x[i]);
481 static void write_em_traj(FILE *fplog, t_commrec *cr,
482 gmx_mdoutf_t outf,
483 gmx_bool bX, gmx_bool bF, const char *confout,
484 gmx_mtop_t *top_global,
485 t_inputrec *ir, gmx_int64_t step,
486 em_state_t *state,
487 t_state *state_global, rvec *f_global)
489 int mdof_flags;
490 gmx_bool bIMDout = FALSE;
493 /* Shall we do IMD output? */
494 if (ir->bIMD)
496 bIMDout = do_per_step(step, IMD_get_step(ir->imd->setup));
499 if ((bX || bF || bIMDout || confout != NULL) && !DOMAINDECOMP(cr))
501 copy_em_coords(state, state_global);
502 f_global = state->f;
505 mdof_flags = 0;
506 if (bX)
508 mdof_flags |= MDOF_X;
510 if (bF)
512 mdof_flags |= MDOF_F;
515 /* If we want IMD output, set appropriate MDOF flag */
516 if (ir->bIMD)
518 mdof_flags |= MDOF_IMD;
521 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
522 top_global, step, (double)step,
523 &state->s, state_global, state->f, f_global);
525 if (confout != NULL && MASTER(cr))
527 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
529 /* Make molecules whole only for confout writing */
530 do_pbc_mtop(fplog, ir->ePBC, state_global->box, top_global,
531 state_global->x);
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,
541 gmx_bool bMolPBC,
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,
545 gmx_int64_t count)
548 t_state *s1, *s2;
549 int i;
550 int start, end;
551 rvec *x1, *x2;
552 real dvdl_constr;
554 s1 = &ems1->s;
555 s2 = &ems2->s;
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);
584 start = 0;
585 end = md->homenr;
587 x1 = s1->x;
588 x2 = s2->x;
590 #pragma omp parallel num_threads(gmx_omp_nthreads_get(emntUpdate))
592 int gf, i, m;
594 gf = 0;
595 #pragma omp for schedule(static) nowait
596 for (i = start; i < end; i++)
598 if (md->cFREEZE)
600 gf = md->cFREEZE[i];
602 for (m = 0; m < DIM; m++)
604 if (ir->opts.nFreeze[gf][m])
606 x2[i][m] = x1[i][m];
608 else
610 x2[i][m] = x1[i][m] + a*f[i][m];
615 if (s2->flags & (1<<estCGP))
617 /* Copy the CG p vector */
618 x1 = s1->cg_p;
619 x2 = s2->cg_p;
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)
632 #pragma omp barrier
633 s2->cg_gl_nalloc = s1->cg_gl_nalloc;
634 srenew(s2->cg_gl, s2->cg_gl_nalloc);
635 #pragma omp barrier
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;
647 if (constr)
649 wallcycle_start(wcycle, ewcCONSTR);
650 dvdl_constr = 0;
651 constrain(NULL, TRUE, TRUE, constr, &top->idef,
652 ir, NULL, cr, count, 0, 1.0, md,
653 s1->x, s2->x, NULL, bMolPBC, s2->box,
654 s2->lambda[efptBONDED], &dvdl_constr,
655 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
656 wallcycle_stop(wcycle, ewcCONSTR);
660 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
661 gmx_mtop_t *top_global, t_inputrec *ir,
662 em_state_t *ems, gmx_localtop_t *top,
663 t_mdatoms *mdatoms, t_forcerec *fr,
664 gmx_vsite_t *vsite, gmx_constr_t constr,
665 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
667 /* Repartition the domain decomposition */
668 wallcycle_start(wcycle, ewcDOMDEC);
669 dd_partition_system(fplog, step, cr, FALSE, 1,
670 NULL, top_global, ir,
671 &ems->s, &ems->f,
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,
685 t_fcdata *fcd,
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)
691 real t;
692 gmx_bool bNS;
693 int nabnsb;
694 tensor force_vir, shake_vir, ekin;
695 real dvdl_constr, prescorr, enercorr, dvdlcorr;
696 real terminate = 0;
698 /* Set the time to the initial time, the time does not change during EM */
699 t = inputrec->init_t;
701 if (bFirst ||
702 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
704 /* This is the first state or an old state used before the last ns */
705 bNS = TRUE;
707 else
709 bNS = FALSE;
710 if (inputrec->nstlist > 0)
712 bNS = TRUE;
714 else if (inputrec->nstlist == -1)
716 nabnsb = natoms_beyond_ns_buffer(inputrec, fr, &top->cgs, NULL, ems->s.x);
717 if (PAR(cr))
719 gmx_sumi(1, &nabnsb, cr);
721 bNS = (nabnsb > 0);
725 if (vsite)
727 construct_vsites(vsite, ems->s.x, 1, NULL,
728 top->idef.iparams, top->idef.il,
729 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
732 if (DOMAINDECOMP(cr) && bNS)
734 /* Repartition the domain decomposition */
735 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
736 ems, top, mdatoms, fr, vsite, constr,
737 nrnb, wcycle);
740 /* Calc force & energy on new trial position */
741 /* do_force always puts the charge groups in the box and shifts again
742 * We do not unshift, so molecules are always whole in congrad.c
744 do_force(fplog, cr, inputrec,
745 count, nrnb, wcycle, top, &top_global->groups,
746 ems->s.box, ems->s.x, &ems->s.hist,
747 ems->f, force_vir, mdatoms, enerd, fcd,
748 ems->s.lambda, graph, fr, vsite, mu_tot, t, NULL, NULL, TRUE,
749 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
750 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
751 (bNS ? GMX_FORCE_NS | GMX_FORCE_DO_LR : 0));
753 /* Clear the unused shake virial and pressure */
754 clear_mat(shake_vir);
755 clear_mat(pres);
757 /* Communicate stuff when parallel */
758 if (PAR(cr) && inputrec->eI != eiNM)
760 wallcycle_start(wcycle, ewcMoveE);
762 global_stat(fplog, gstat, cr, enerd, force_vir, shake_vir, mu_tot,
763 inputrec, NULL, NULL, NULL, 1, &terminate,
764 top_global, &ems->s, FALSE,
765 CGLO_ENERGY |
766 CGLO_PRESSURE |
767 CGLO_CONSTRAINT |
768 CGLO_FIRSTITERATE);
770 wallcycle_stop(wcycle, ewcMoveE);
773 /* Calculate long range corrections to pressure and energy */
774 calc_dispcorr(fplog, inputrec, fr, count, top_global->natoms, ems->s.box, ems->s.lambda[efptVDW],
775 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
776 enerd->term[F_DISPCORR] = enercorr;
777 enerd->term[F_EPOT] += enercorr;
778 enerd->term[F_PRES] += prescorr;
779 enerd->term[F_DVDL] += dvdlcorr;
781 ems->epot = enerd->term[F_EPOT];
783 if (constr)
785 /* Project out the constraint components of the force */
786 wallcycle_start(wcycle, ewcCONSTR);
787 dvdl_constr = 0;
788 constrain(NULL, FALSE, FALSE, constr, &top->idef,
789 inputrec, NULL, cr, count, 0, 1.0, mdatoms,
790 ems->s.x, ems->f, ems->f, fr->bMolPBC, ems->s.box,
791 ems->s.lambda[efptBONDED], &dvdl_constr,
792 NULL, &shake_vir, nrnb, econqForceDispl, FALSE, 0, 0);
793 if (fr->bSepDVDL && fplog)
795 gmx_print_sepdvdl(fplog, "Constraints", t, dvdl_constr);
797 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
798 m_add(force_vir, shake_vir, vir);
799 wallcycle_stop(wcycle, ewcCONSTR);
801 else
803 copy_mat(force_vir, vir);
806 clear_mat(ekin);
807 enerd->term[F_PRES] =
808 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
810 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
812 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
814 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
818 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
819 gmx_mtop_t *mtop,
820 em_state_t *s_min, em_state_t *s_b)
822 rvec *fm, *fb, *fmg;
823 t_block *cgs_gl;
824 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
825 double partsum;
826 unsigned char *grpnrFREEZE;
828 if (debug)
830 fprintf(debug, "Doing reorder_partsum\n");
833 fm = s_min->f;
834 fb = s_b->f;
836 cgs_gl = dd_charge_groups_global(cr->dd);
837 index = cgs_gl->index;
839 /* Collect fm in a global vector fmg.
840 * This conflicts with the spirit of domain decomposition,
841 * but to fully optimize this a much more complicated algorithm is required.
843 snew(fmg, mtop->natoms);
845 ncg = s_min->s.ncg_gl;
846 cg_gl = s_min->s.cg_gl;
847 i = 0;
848 for (c = 0; c < ncg; c++)
850 cg = cg_gl[c];
851 a0 = index[cg];
852 a1 = index[cg+1];
853 for (a = a0; a < a1; a++)
855 copy_rvec(fm[i], fmg[a]);
856 i++;
859 gmx_sum(mtop->natoms*3, fmg[0], cr);
861 /* Now we will determine the part of the sum for the cgs in state s_b */
862 ncg = s_b->s.ncg_gl;
863 cg_gl = s_b->s.cg_gl;
864 partsum = 0;
865 i = 0;
866 gf = 0;
867 grpnrFREEZE = mtop->groups.grpnr[egcFREEZE];
868 for (c = 0; c < ncg; c++)
870 cg = cg_gl[c];
871 a0 = index[cg];
872 a1 = index[cg+1];
873 for (a = a0; a < a1; a++)
875 if (mdatoms->cFREEZE && grpnrFREEZE)
877 gf = grpnrFREEZE[i];
879 for (m = 0; m < DIM; m++)
881 if (!opts->nFreeze[gf][m])
883 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
886 i++;
890 sfree(fmg);
892 return partsum;
895 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
896 gmx_mtop_t *mtop,
897 em_state_t *s_min, em_state_t *s_b)
899 rvec *fm, *fb;
900 double sum;
901 int gf, i, m;
903 /* This is just the classical Polak-Ribiere calculation of beta;
904 * it looks a bit complicated since we take freeze groups into account,
905 * and might have to sum it in parallel runs.
908 if (!DOMAINDECOMP(cr) ||
909 (s_min->s.ddp_count == cr->dd->ddp_count &&
910 s_b->s.ddp_count == cr->dd->ddp_count))
912 fm = s_min->f;
913 fb = s_b->f;
914 sum = 0;
915 gf = 0;
916 /* This part of code can be incorrect with DD,
917 * since the atom ordering in s_b and s_min might differ.
919 for (i = 0; i < mdatoms->homenr; i++)
921 if (mdatoms->cFREEZE)
923 gf = mdatoms->cFREEZE[i];
925 for (m = 0; m < DIM; m++)
927 if (!opts->nFreeze[gf][m])
929 sum += (fb[i][m] - fm[i][m])*fb[i][m];
934 else
936 /* We need to reorder cgs while summing */
937 sum = reorder_partsum(cr, opts, mdatoms, mtop, s_min, s_b);
939 if (PAR(cr))
941 gmx_sumd(1, &sum, cr);
944 return sum/sqr(s_min->fnorm);
947 double do_cg(FILE *fplog, t_commrec *cr,
948 int nfile, const t_filenm fnm[],
949 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
950 int gmx_unused nstglobalcomm,
951 gmx_vsite_t *vsite, gmx_constr_t constr,
952 int gmx_unused stepout,
953 t_inputrec *inputrec,
954 gmx_mtop_t *top_global, t_fcdata *fcd,
955 t_state *state_global,
956 t_mdatoms *mdatoms,
957 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
958 gmx_edsam_t gmx_unused ed,
959 t_forcerec *fr,
960 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
961 gmx_membed_t gmx_unused membed,
962 real gmx_unused cpt_period, real gmx_unused max_hours,
963 const char gmx_unused *deviceOptions,
964 int imdport,
965 unsigned long gmx_unused Flags,
966 gmx_walltime_accounting_t walltime_accounting)
968 const char *CG = "Polak-Ribiere Conjugate Gradients";
970 em_state_t *s_min, *s_a, *s_b, *s_c;
971 gmx_localtop_t *top;
972 gmx_enerdata_t *enerd;
973 rvec *f;
974 gmx_global_stat_t gstat;
975 t_graph *graph;
976 rvec *f_global, *p, *sf, *sfm;
977 double gpa, gpb, gpc, tmp, sum[2], minstep;
978 real fnormn;
979 real stepsize;
980 real a, b, c, beta = 0.0;
981 real epot_repl = 0;
982 real pnorm;
983 t_mdebin *mdebin;
984 gmx_bool converged, foundlower;
985 rvec mu_tot;
986 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
987 tensor vir, pres;
988 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
989 gmx_mdoutf_t outf;
990 int i, m, gf, step, nminstep;
991 real terminate = 0;
993 step = 0;
995 s_min = init_em_state();
996 s_a = init_em_state();
997 s_b = init_em_state();
998 s_c = init_em_state();
1000 /* Init em and store the local state in s_min */
1001 init_em(fplog, CG, cr, inputrec,
1002 state_global, top_global, s_min, &top, &f, &f_global,
1003 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1004 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
1006 /* Print to log file */
1007 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1009 /* Max number of steps */
1010 number_steps = inputrec->nsteps;
1012 if (MASTER(cr))
1014 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1016 if (fplog)
1018 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1021 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1022 /* do_force always puts the charge groups in the box and shifts again
1023 * We do not unshift, so molecules are always whole in congrad.c
1025 evaluate_energy(fplog, cr,
1026 top_global, s_min, top,
1027 inputrec, nrnb, wcycle, gstat,
1028 vsite, constr, fcd, graph, mdatoms, fr,
1029 mu_tot, enerd, vir, pres, -1, TRUE);
1030 where();
1032 if (MASTER(cr))
1034 /* Copy stuff to the energy bin for easy printing etc. */
1035 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1036 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1037 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1039 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1040 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1041 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1043 where();
1045 /* Estimate/guess the initial stepsize */
1046 stepsize = inputrec->em_stepsize/s_min->fnorm;
1048 if (MASTER(cr))
1050 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1051 s_min->fmax, s_min->a_fmax+1);
1052 fprintf(stderr, " F-Norm = %12.5e\n",
1053 s_min->fnorm/sqrt(state_global->natoms));
1054 fprintf(stderr, "\n");
1055 /* and copy to the log file too... */
1056 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1057 s_min->fmax, s_min->a_fmax+1);
1058 fprintf(fplog, " F-Norm = %12.5e\n",
1059 s_min->fnorm/sqrt(state_global->natoms));
1060 fprintf(fplog, "\n");
1062 /* Start the loop over CG steps.
1063 * Each successful step is counted, and we continue until
1064 * we either converge or reach the max number of steps.
1066 converged = FALSE;
1067 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1070 /* start taking steps in a new direction
1071 * First time we enter the routine, beta=0, and the direction is
1072 * simply the negative gradient.
1075 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1076 p = s_min->s.cg_p;
1077 sf = s_min->f;
1078 gpa = 0;
1079 gf = 0;
1080 for (i = 0; i < mdatoms->homenr; i++)
1082 if (mdatoms->cFREEZE)
1084 gf = mdatoms->cFREEZE[i];
1086 for (m = 0; m < DIM; m++)
1088 if (!inputrec->opts.nFreeze[gf][m])
1090 p[i][m] = sf[i][m] + beta*p[i][m];
1091 gpa -= p[i][m]*sf[i][m];
1092 /* f is negative gradient, thus the sign */
1094 else
1096 p[i][m] = 0;
1101 /* Sum the gradient along the line across CPUs */
1102 if (PAR(cr))
1104 gmx_sumd(1, &gpa, cr);
1107 /* Calculate the norm of the search vector */
1108 get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1110 /* Just in case stepsize reaches zero due to numerical precision... */
1111 if (stepsize <= 0)
1113 stepsize = inputrec->em_stepsize/pnorm;
1117 * Double check the value of the derivative in the search direction.
1118 * If it is positive it must be due to the old information in the
1119 * CG formula, so just remove that and start over with beta=0.
1120 * This corresponds to a steepest descent step.
1122 if (gpa > 0)
1124 beta = 0;
1125 step--; /* Don't count this step since we are restarting */
1126 continue; /* Go back to the beginning of the big for-loop */
1129 /* Calculate minimum allowed stepsize, before the average (norm)
1130 * relative change in coordinate is smaller than precision
1132 minstep = 0;
1133 for (i = 0; i < mdatoms->homenr; i++)
1135 for (m = 0; m < DIM; m++)
1137 tmp = fabs(s_min->s.x[i][m]);
1138 if (tmp < 1.0)
1140 tmp = 1.0;
1142 tmp = p[i][m]/tmp;
1143 minstep += tmp*tmp;
1146 /* Add up from all CPUs */
1147 if (PAR(cr))
1149 gmx_sumd(1, &minstep, cr);
1152 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1154 if (stepsize < minstep)
1156 converged = TRUE;
1157 break;
1160 /* Write coordinates if necessary */
1161 do_x = do_per_step(step, inputrec->nstxout);
1162 do_f = do_per_step(step, inputrec->nstfout);
1164 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1165 top_global, inputrec, step,
1166 s_min, state_global, f_global);
1168 /* Take a step downhill.
1169 * In theory, we should minimize the function along this direction.
1170 * That is quite possible, but it turns out to take 5-10 function evaluations
1171 * for each line. However, we dont really need to find the exact minimum -
1172 * it is much better to start a new CG step in a modified direction as soon
1173 * as we are close to it. This will save a lot of energy evaluations.
1175 * In practice, we just try to take a single step.
1176 * If it worked (i.e. lowered the energy), we increase the stepsize but
1177 * the continue straight to the next CG step without trying to find any minimum.
1178 * If it didn't work (higher energy), there must be a minimum somewhere between
1179 * the old position and the new one.
1181 * Due to the finite numerical accuracy, it turns out that it is a good idea
1182 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1183 * This leads to lower final energies in the tests I've done. / Erik
1185 s_a->epot = s_min->epot;
1186 a = 0.0;
1187 c = a + stepsize; /* reference position along line is zero */
1189 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1191 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1192 s_min, top, mdatoms, fr, vsite, constr,
1193 nrnb, wcycle);
1196 /* Take a trial step (new coords in s_c) */
1197 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1198 constr, top, nrnb, wcycle, -1);
1200 neval++;
1201 /* Calculate energy for the trial step */
1202 evaluate_energy(fplog, cr,
1203 top_global, s_c, top,
1204 inputrec, nrnb, wcycle, gstat,
1205 vsite, constr, fcd, graph, mdatoms, fr,
1206 mu_tot, enerd, vir, pres, -1, FALSE);
1208 /* Calc derivative along line */
1209 p = s_c->s.cg_p;
1210 sf = s_c->f;
1211 gpc = 0;
1212 for (i = 0; i < mdatoms->homenr; i++)
1214 for (m = 0; m < DIM; m++)
1216 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1219 /* Sum the gradient along the line across CPUs */
1220 if (PAR(cr))
1222 gmx_sumd(1, &gpc, cr);
1225 /* This is the max amount of increase in energy we tolerate */
1226 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1228 /* Accept the step if the energy is lower, or if it is not significantly higher
1229 * and the line derivative is still negative.
1231 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1233 foundlower = TRUE;
1234 /* Great, we found a better energy. Increase step for next iteration
1235 * if we are still going down, decrease it otherwise
1237 if (gpc < 0)
1239 stepsize *= 1.618034; /* The golden section */
1241 else
1243 stepsize *= 0.618034; /* 1/golden section */
1246 else
1248 /* New energy is the same or higher. We will have to do some work
1249 * to find a smaller value in the interval. Take smaller step next time!
1251 foundlower = FALSE;
1252 stepsize *= 0.618034;
1258 /* OK, if we didn't find a lower value we will have to locate one now - there must
1259 * be one in the interval [a=0,c].
1260 * The same thing is valid here, though: Don't spend dozens of iterations to find
1261 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1262 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1264 * I also have a safeguard for potentially really patological functions so we never
1265 * take more than 20 steps before we give up ...
1267 * If we already found a lower value we just skip this step and continue to the update.
1269 if (!foundlower)
1271 nminstep = 0;
1275 /* Select a new trial point.
1276 * If the derivatives at points a & c have different sign we interpolate to zero,
1277 * otherwise just do a bisection.
1279 if (gpa < 0 && gpc > 0)
1281 b = a + gpa*(a-c)/(gpc-gpa);
1283 else
1285 b = 0.5*(a+c);
1288 /* safeguard if interpolation close to machine accuracy causes errors:
1289 * never go outside the interval
1291 if (b <= a || b >= c)
1293 b = 0.5*(a+c);
1296 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1298 /* Reload the old state */
1299 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1300 s_min, top, mdatoms, fr, vsite, constr,
1301 nrnb, wcycle);
1304 /* Take a trial step to this new point - new coords in s_b */
1305 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1306 constr, top, nrnb, wcycle, -1);
1308 neval++;
1309 /* Calculate energy for the trial step */
1310 evaluate_energy(fplog, cr,
1311 top_global, s_b, top,
1312 inputrec, nrnb, wcycle, gstat,
1313 vsite, constr, fcd, graph, mdatoms, fr,
1314 mu_tot, enerd, vir, pres, -1, FALSE);
1316 /* p does not change within a step, but since the domain decomposition
1317 * might change, we have to use cg_p of s_b here.
1319 p = s_b->s.cg_p;
1320 sf = s_b->f;
1321 gpb = 0;
1322 for (i = 0; i < mdatoms->homenr; i++)
1324 for (m = 0; m < DIM; m++)
1326 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1329 /* Sum the gradient along the line across CPUs */
1330 if (PAR(cr))
1332 gmx_sumd(1, &gpb, cr);
1335 if (debug)
1337 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1338 s_a->epot, s_b->epot, s_c->epot, gpb);
1341 epot_repl = s_b->epot;
1343 /* Keep one of the intervals based on the value of the derivative at the new point */
1344 if (gpb > 0)
1346 /* Replace c endpoint with b */
1347 swap_em_state(s_b, s_c);
1348 c = b;
1349 gpc = gpb;
1351 else
1353 /* Replace a endpoint with b */
1354 swap_em_state(s_b, s_a);
1355 a = b;
1356 gpa = gpb;
1360 * Stop search as soon as we find a value smaller than the endpoints.
1361 * Never run more than 20 steps, no matter what.
1363 nminstep++;
1365 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1366 (nminstep < 20));
1368 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1369 nminstep >= 20)
1371 /* OK. We couldn't find a significantly lower energy.
1372 * If beta==0 this was steepest descent, and then we give up.
1373 * If not, set beta=0 and restart with steepest descent before quitting.
1375 if (beta == 0.0)
1377 /* Converged */
1378 converged = TRUE;
1379 break;
1381 else
1383 /* Reset memory before giving up */
1384 beta = 0.0;
1385 continue;
1389 /* Select min energy state of A & C, put the best in B.
1391 if (s_c->epot < s_a->epot)
1393 if (debug)
1395 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1396 s_c->epot, s_a->epot);
1398 swap_em_state(s_b, s_c);
1399 gpb = gpc;
1400 b = c;
1402 else
1404 if (debug)
1406 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1407 s_a->epot, s_c->epot);
1409 swap_em_state(s_b, s_a);
1410 gpb = gpa;
1411 b = a;
1415 else
1417 if (debug)
1419 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1420 s_c->epot);
1422 swap_em_state(s_b, s_c);
1423 gpb = gpc;
1424 b = c;
1427 /* new search direction */
1428 /* beta = 0 means forget all memory and restart with steepest descents. */
1429 if (nstcg && ((step % nstcg) == 0))
1431 beta = 0.0;
1433 else
1435 /* s_min->fnorm cannot be zero, because then we would have converged
1436 * and broken out.
1439 /* Polak-Ribiere update.
1440 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1442 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1444 /* Limit beta to prevent oscillations */
1445 if (fabs(beta) > 5.0)
1447 beta = 0.0;
1451 /* update positions */
1452 swap_em_state(s_min, s_b);
1453 gpa = gpb;
1455 /* Print it if necessary */
1456 if (MASTER(cr))
1458 if (bVerbose)
1460 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1461 step, s_min->epot, s_min->fnorm/sqrt(state_global->natoms),
1462 s_min->fmax, s_min->a_fmax+1);
1464 /* Store the new (lower) energies */
1465 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1466 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1467 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1469 do_log = do_per_step(step, inputrec->nstlog);
1470 do_ene = do_per_step(step, inputrec->nstenergy);
1472 /* Prepare IMD energy record, if bIMD is TRUE. */
1473 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1475 if (do_log)
1477 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1479 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1480 do_log ? fplog : NULL, step, step, eprNORMAL,
1481 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1484 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1485 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
1487 IMD_send_positions(inputrec->imd);
1490 /* Stop when the maximum force lies below tolerance.
1491 * If we have reached machine precision, converged is already set to true.
1493 converged = converged || (s_min->fmax < inputrec->em_tol);
1495 } /* End of the loop */
1497 /* IMD cleanup, if bIMD is TRUE. */
1498 IMD_finalize(inputrec->bIMD, inputrec->imd);
1500 if (converged)
1502 step--; /* we never took that last step in this case */
1505 if (s_min->fmax > inputrec->em_tol)
1507 if (MASTER(cr))
1509 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1510 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1512 converged = FALSE;
1515 if (MASTER(cr))
1517 /* If we printed energy and/or logfile last step (which was the last step)
1518 * we don't have to do it again, but otherwise print the final values.
1520 if (!do_log)
1522 /* Write final value to log since we didn't do anything the last step */
1523 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1525 if (!do_ene || !do_log)
1527 /* Write final energy file entries */
1528 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1529 !do_log ? fplog : NULL, step, step, eprNORMAL,
1530 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1534 /* Print some stuff... */
1535 if (MASTER(cr))
1537 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1540 /* IMPORTANT!
1541 * For accurate normal mode calculation it is imperative that we
1542 * store the last conformation into the full precision binary trajectory.
1544 * However, we should only do it if we did NOT already write this step
1545 * above (which we did if do_x or do_f was true).
1547 do_x = !do_per_step(step, inputrec->nstxout);
1548 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1550 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1551 top_global, inputrec, step,
1552 s_min, state_global, f_global);
1554 fnormn = s_min->fnorm/sqrt(state_global->natoms);
1556 if (MASTER(cr))
1558 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1559 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1560 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1561 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1563 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1566 finish_em(cr, outf, walltime_accounting, wcycle);
1568 /* To print the actual number of steps we needed somewhere */
1569 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1571 return 0;
1572 } /* That's all folks */
1575 double do_lbfgs(FILE *fplog, t_commrec *cr,
1576 int nfile, const t_filenm fnm[],
1577 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
1578 int gmx_unused nstglobalcomm,
1579 gmx_vsite_t *vsite, gmx_constr_t constr,
1580 int gmx_unused stepout,
1581 t_inputrec *inputrec,
1582 gmx_mtop_t *top_global, t_fcdata *fcd,
1583 t_state *state,
1584 t_mdatoms *mdatoms,
1585 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1586 gmx_edsam_t gmx_unused ed,
1587 t_forcerec *fr,
1588 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1589 gmx_membed_t gmx_unused membed,
1590 real gmx_unused cpt_period, real gmx_unused max_hours,
1591 const char gmx_unused *deviceOptions,
1592 int imdport,
1593 unsigned long gmx_unused Flags,
1594 gmx_walltime_accounting_t walltime_accounting)
1596 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1597 em_state_t ems;
1598 gmx_localtop_t *top;
1599 gmx_enerdata_t *enerd;
1600 rvec *f;
1601 gmx_global_stat_t gstat;
1602 t_graph *graph;
1603 rvec *f_global;
1604 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1605 double stepsize, gpa, gpb, gpc, tmp, minstep;
1606 real *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1607 real *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1608 real a, b, c, maxdelta, delta;
1609 real diag, Epot0, Epot, EpotA, EpotB, EpotC;
1610 real dgdx, dgdg, sq, yr, beta;
1611 t_mdebin *mdebin;
1612 gmx_bool converged, first;
1613 rvec mu_tot;
1614 real fnorm, fmax;
1615 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1616 tensor vir, pres;
1617 int start, end, number_steps;
1618 gmx_mdoutf_t outf;
1619 int i, k, m, n, nfmax, gf, step;
1620 int mdof_flags;
1621 /* not used */
1622 real terminate;
1624 if (PAR(cr))
1626 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1629 if (NULL != constr)
1631 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).");
1634 n = 3*state->natoms;
1635 nmaxcorr = inputrec->nbfgscorr;
1637 /* Allocate memory */
1638 /* Use pointers to real so we dont have to loop over both atoms and
1639 * dimensions all the time...
1640 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1641 * that point to the same memory.
1643 snew(xa, n);
1644 snew(xb, n);
1645 snew(xc, n);
1646 snew(fa, n);
1647 snew(fb, n);
1648 snew(fc, n);
1649 snew(frozen, n);
1651 snew(p, n);
1652 snew(lastx, n);
1653 snew(lastf, n);
1654 snew(rho, nmaxcorr);
1655 snew(alpha, nmaxcorr);
1657 snew(dx, nmaxcorr);
1658 for (i = 0; i < nmaxcorr; i++)
1660 snew(dx[i], n);
1663 snew(dg, nmaxcorr);
1664 for (i = 0; i < nmaxcorr; i++)
1666 snew(dg[i], n);
1669 step = 0;
1670 neval = 0;
1672 /* Init em */
1673 init_em(fplog, LBFGS, cr, inputrec,
1674 state, top_global, &ems, &top, &f, &f_global,
1675 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1676 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
1677 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1678 * so we free some memory again.
1680 sfree(ems.s.x);
1681 sfree(ems.f);
1683 xx = (real *)state->x;
1684 ff = (real *)f;
1686 start = 0;
1687 end = mdatoms->homenr;
1689 /* Print to log file */
1690 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1692 do_log = do_ene = do_x = do_f = TRUE;
1694 /* Max number of steps */
1695 number_steps = inputrec->nsteps;
1697 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1698 gf = 0;
1699 for (i = start; i < end; i++)
1701 if (mdatoms->cFREEZE)
1703 gf = mdatoms->cFREEZE[i];
1705 for (m = 0; m < DIM; m++)
1707 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1710 if (MASTER(cr))
1712 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1714 if (fplog)
1716 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1719 if (vsite)
1721 construct_vsites(vsite, state->x, 1, NULL,
1722 top->idef.iparams, top->idef.il,
1723 fr->ePBC, fr->bMolPBC, cr, state->box);
1726 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1727 /* do_force always puts the charge groups in the box and shifts again
1728 * We do not unshift, so molecules are always whole
1730 neval++;
1731 ems.s.x = state->x;
1732 ems.f = f;
1733 evaluate_energy(fplog, cr,
1734 top_global, &ems, top,
1735 inputrec, nrnb, wcycle, gstat,
1736 vsite, constr, fcd, graph, mdatoms, fr,
1737 mu_tot, enerd, vir, pres, -1, TRUE);
1738 where();
1740 if (MASTER(cr))
1742 /* Copy stuff to the energy bin for easy printing etc. */
1743 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1744 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
1745 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1747 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
1748 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1749 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1751 where();
1753 /* This is the starting energy */
1754 Epot = enerd->term[F_EPOT];
1756 fnorm = ems.fnorm;
1757 fmax = ems.fmax;
1758 nfmax = ems.a_fmax;
1760 /* Set the initial step.
1761 * since it will be multiplied by the non-normalized search direction
1762 * vector (force vector the first time), we scale it by the
1763 * norm of the force.
1766 if (MASTER(cr))
1768 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1769 fprintf(stderr, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1770 fprintf(stderr, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1771 fprintf(stderr, "\n");
1772 /* and copy to the log file too... */
1773 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1774 fprintf(fplog, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1775 fprintf(fplog, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1776 fprintf(fplog, "\n");
1779 point = 0;
1780 for (i = 0; i < n; i++)
1782 if (!frozen[i])
1784 dx[point][i] = ff[i]; /* Initial search direction */
1786 else
1788 dx[point][i] = 0;
1792 stepsize = 1.0/fnorm;
1793 converged = FALSE;
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.
1800 ncorr = 0;
1802 /* Set the gradient from the force */
1803 converged = FALSE;
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);
1811 mdof_flags = 0;
1812 if (do_x)
1814 mdof_flags |= MDOF_X;
1817 if (do_f)
1819 mdof_flags |= MDOF_F;
1822 if (inputrec->bIMD)
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 */
1833 s = dx[point];
1835 /* calculate line gradient */
1836 for (gpa = 0, i = 0; i < n; i++)
1838 gpa -= s[i]*ff[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++)
1846 tmp = fabs(xx[i]);
1847 if (tmp < 1.0)
1849 tmp = 1.0;
1851 tmp = s[i]/tmp;
1852 minstep += tmp*tmp;
1854 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1856 if (stepsize < minstep)
1858 converged = TRUE;
1859 break;
1862 /* Store old forces and coordinates */
1863 for (i = 0; i < n; i++)
1865 lastx[i] = xx[i];
1866 lastf[i] = ff[i];
1868 Epot0 = Epot;
1870 first = TRUE;
1872 for (i = 0; i < n; i++)
1874 xa[i] = xx[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
1894 foundlower = FALSE;
1895 EpotA = Epot0;
1896 a = 0.0;
1897 c = a + stepsize; /* reference position along line is zero */
1899 /* Check stepsize first. We do not allow displacements
1900 * larger than emstep.
1904 c = a + stepsize;
1905 maxdelta = 0;
1906 for (i = 0; i < n; i++)
1908 delta = c*s[i];
1909 if (delta > maxdelta)
1911 maxdelta = delta;
1914 if (maxdelta > inputrec->em_stepsize)
1916 stepsize *= 0.1;
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];
1927 neval++;
1928 /* Calculate energy for the trial step */
1929 ems.s.x = (rvec *)xc;
1930 ems.f = (rvec *)fc;
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);
1936 EpotC = ems.epot;
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 */
1944 if (PAR(cr))
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)))
1957 foundlower = TRUE;
1958 /* Great, we found a better energy. Increase step for next iteration
1959 * if we are still going down, decrease it otherwise
1961 if (gpc < 0)
1963 stepsize *= 1.618034; /* The golden section */
1965 else
1967 stepsize *= 0.618034; /* 1/golden section */
1970 else
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!
1975 foundlower = FALSE;
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.
1991 if (!foundlower)
1994 nminstep = 0;
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);
2006 else
2008 b = 0.5*(a+c);
2011 /* safeguard if interpolation close to machine accuracy causes errors:
2012 * never go outside the interval
2014 if (b <= a || b >= c)
2016 b = 0.5*(a+c);
2019 /* Take a trial step */
2020 for (i = 0; i < n; i++)
2022 xb[i] = lastx[i] + b*s[i];
2025 neval++;
2026 /* Calculate energy for the trial step */
2027 ems.s.x = (rvec *)xb;
2028 ems.f = (rvec *)fb;
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);
2034 EpotB = ems.epot;
2036 fnorm = ems.fnorm;
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 */
2044 if (PAR(cr))
2046 gmx_sumd(1, &gpb, cr);
2049 /* Keep one of the intervals based on the value of the derivative at the new point */
2050 if (gpb > 0)
2052 /* Replace c endpoint with b */
2053 EpotC = EpotB;
2054 c = b;
2055 gpc = gpb;
2056 /* swap coord pointers b/c */
2057 xtmp = xb;
2058 ftmp = fb;
2059 xb = xc;
2060 fb = fc;
2061 xc = xtmp;
2062 fc = ftmp;
2064 else
2066 /* Replace a endpoint with b */
2067 EpotA = EpotB;
2068 a = b;
2069 gpa = gpb;
2070 /* swap coord pointers a/b */
2071 xtmp = xb;
2072 ftmp = fb;
2073 xb = xa;
2074 fb = fa;
2075 xa = xtmp;
2076 fa = ftmp;
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.
2084 nminstep++;
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.
2094 if (ncorr == 0)
2096 /* Converged */
2097 converged = TRUE;
2098 break;
2100 else
2102 /* Reset memory */
2103 ncorr = 0;
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;
2111 continue;
2115 /* Select min energy state of A & C, put the best in xx/ff/Epot
2117 if (EpotC < EpotA)
2119 Epot = EpotC;
2120 /* Use state C */
2121 for (i = 0; i < n; i++)
2123 xx[i] = xc[i];
2124 ff[i] = fc[i];
2126 stepsize = c;
2128 else
2130 Epot = EpotA;
2131 /* Use state A */
2132 for (i = 0; i < n; i++)
2134 xx[i] = xa[i];
2135 ff[i] = fa[i];
2137 stepsize = a;
2141 else
2143 /* found lower */
2144 Epot = EpotC;
2145 /* Use state C */
2146 for (i = 0; i < n; i++)
2148 xx[i] = xc[i];
2149 ff[i] = fc[i];
2151 stepsize = c;
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)
2161 ncorr++;
2164 for (i = 0; i < n; i++)
2166 dg[point][i] = lastf[i]-ff[i];
2167 dx[point][i] *= stepsize;
2170 dgdg = 0;
2171 dgdx = 0;
2172 for (i = 0; i < n; i++)
2174 dgdg += dg[point][i]*dg[point][i];
2175 dgdx += dg[point][i]*dx[point][i];
2178 diag = dgdx/dgdg;
2180 rho[point] = 1.0/dgdx;
2181 point++;
2183 if (point >= nmaxcorr)
2185 point = 0;
2188 /* Update */
2189 for (i = 0; i < n; i++)
2191 p[i] = ff[i];
2194 cp = point;
2196 /* Recursive update. First go back over the memory points */
2197 for (k = 0; k < ncorr; k++)
2199 cp--;
2200 if (cp < 0)
2202 cp = ncorr-1;
2205 sq = 0;
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++)
2221 p[i] *= diag;
2224 /* And then go forward again */
2225 for (k = 0; k < ncorr; k++)
2227 yr = 0;
2228 for (i = 0; i < n; i++)
2230 yr += p[i]*dg[cp][i];
2233 beta = rho[cp]*yr;
2234 beta = alpha[cp]-beta;
2236 for (i = 0; i < n; i++)
2238 p[i] += beta*dx[cp][i];
2241 cp++;
2242 if (cp >= ncorr)
2244 cp = 0;
2248 for (i = 0; i < n; i++)
2250 if (!frozen[i])
2252 dx[point][i] = p[i];
2254 else
2256 dx[point][i] = 0;
2260 stepsize = 1.0;
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 */
2266 if (MASTER(cr))
2268 if (bVerbose)
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);
2279 if (do_log)
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);
2305 if (converged)
2307 step--; /* we never took that last step in this case */
2310 if (fmax > inputrec->em_tol)
2312 if (MASTER(cr))
2314 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2315 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2317 converged = 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... */
2335 if (MASTER(cr))
2337 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2340 /* IMPORTANT!
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,
2351 &ems, state, f);
2353 if (MASTER(cr))
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);
2368 return 0;
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,
2381 t_mdatoms *mdatoms,
2382 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2383 gmx_edsam_t gmx_unused ed,
2384 t_forcerec *fr,
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,
2389 int imdport,
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;
2395 rvec *f_global;
2396 gmx_localtop_t *top;
2397 gmx_enerdata_t *enerd;
2398 rvec *f;
2399 gmx_global_stat_t gstat;
2400 t_graph *graph;
2401 real stepsize, constepsize;
2402 real ustep, fnormn;
2403 gmx_mdoutf_t outf;
2404 t_mdebin *mdebin;
2405 gmx_bool bDone, bAbort, do_x, do_f;
2406 tensor vir, pres;
2407 rvec mu_tot;
2408 int nsteps;
2409 int count = 0;
2410 int steps_accepted = 0;
2411 /* not used */
2412 real terminate = 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, wcycle);
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;
2430 stepsize = 0;
2432 /* Max number of steps */
2433 nsteps = inputrec->nsteps;
2435 if (MASTER(cr))
2437 /* Print to the screen */
2438 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2440 if (fplog)
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
2451 count = 0;
2452 bDone = FALSE;
2453 bAbort = FALSE;
2454 while (!bDone && !bAbort)
2456 bAbort = (nsteps >= 0) && (count == nsteps);
2458 /* set new coordinates, except for first step */
2459 if (count > 0)
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);
2472 if (MASTER(cr))
2474 print_ebin_header(fplog, count, count, s_try->s.lambda[efptFEP]);
2477 if (count == 0)
2479 s_min->epot = s_try->epot + 1;
2482 /* Print it if necessary */
2483 if (MASTER(cr))
2485 if (bVerbose)
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));
2507 fflush(fplog);
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) )
2518 steps_accepted++;
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
2525 sampled energy */
2526 swap_em_state(s_min, s_try);
2527 if (count > 0)
2529 ustep *= 1.2;
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);
2539 else
2541 /* If energy is not smaller make the step smaller... */
2542 ustep *= 0.5;
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,
2549 nrnb, wcycle);
2553 /* Determine new step */
2554 stepsize = ustep/s_min->fmax;
2556 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2557 #ifdef GMX_DOUBLE
2558 if (count == nsteps || ustep < 1e-12)
2559 #else
2560 if (count == nsteps || ustep < 1e-6)
2561 #endif
2563 if (MASTER(cr))
2565 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2566 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2568 bAbort = TRUE;
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);
2577 count++;
2578 } /* End of the loop */
2580 /* IMD cleanup, if bIMD is TRUE. */
2581 IMD_finalize(inputrec->bIMD, inputrec->imd);
2583 /* Print some data... */
2584 if (MASTER(cr))
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);
2594 if (MASTER(cr))
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);
2609 return 0;
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,
2622 t_mdatoms *mdatoms,
2623 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2624 gmx_edsam_t gmx_unused ed,
2625 t_forcerec *fr,
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,
2630 int imdport,
2631 unsigned long gmx_unused Flags,
2632 gmx_walltime_accounting_t walltime_accounting)
2634 const char *NM = "Normal Mode Analysis";
2635 gmx_mdoutf_t outf;
2636 int natoms, atom, d;
2637 int nnodes, node;
2638 rvec *f_global;
2639 gmx_localtop_t *top;
2640 gmx_enerdata_t *enerd;
2641 rvec *f;
2642 gmx_global_stat_t gstat;
2643 t_graph *graph;
2644 real t, t0, lambda, lam0;
2645 gmx_bool bNS;
2646 tensor vir, pres;
2647 rvec mu_tot;
2648 rvec *fneg, *dfdx;
2649 gmx_bool bSparse; /* use sparse matrix storage format */
2650 size_t sz = 0;
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);
2658 real x_min;
2659 real fnorm, fmax;
2661 if (constr != NULL)
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,
2671 &f, &f_global,
2672 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2673 nfile, fnm, &outf, NULL, imdport, Flags, wcycle);
2675 natoms = top_global->natoms;
2676 snew(fneg, natoms);
2677 snew(dfdx, natoms);
2679 #ifndef GMX_DOUBLE
2680 if (MASTER(cr))
2682 fprintf(stderr,
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");
2688 #endif
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");
2699 bSparse = FALSE;
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);
2704 bSparse = FALSE;
2706 else
2708 md_print_info(cr, fplog, "Using compressed symmetric sparse Hessian format.\n");
2709 bSparse = TRUE;
2712 if (MASTER(cr))
2714 sz = DIM*top_global->natoms;
2716 fprintf(stderr, "Allocating Hessian memory...\n\n");
2718 if (bSparse)
2720 sparse_matrix = gmx_sparsematrix_init(sz);
2721 sparse_matrix->compressed_symmetric = TRUE;
2723 else
2725 snew(full_matrix, sz*sz);
2729 /* Initial values */
2730 t0 = inputrec->init_t;
2731 lam0 = inputrec->fepvals->init_lambda;
2732 t = t0;
2733 lambda = lam0;
2735 init_nrnb(nrnb);
2737 where();
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;
2745 if (MASTER(cr))
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 */
2754 cr->nnodes = 1;
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 */
2795 cr->nnodes = 1;
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++)
2823 dfdx[j][k] =
2824 -(state_work->f[j][k] - fneg[j][k])/(2*der_range);
2828 if (!MASTER(cr))
2830 #ifdef GMX_MPI
2831 #ifdef GMX_DOUBLE
2832 #define mpi_type MPI_DOUBLE
2833 #else
2834 #define mpi_type MPI_FLOAT
2835 #endif
2836 MPI_Send(dfdx[0], natoms*DIM, mpi_type, MASTERNODE(cr), cr->nodeid,
2837 cr->mpi_comm_mygroup);
2838 #endif
2840 else
2842 for (node = 0; (node < nnodes && atom+node < natoms); node++)
2844 if (node > 0)
2846 #ifdef GMX_MPI
2847 MPI_Status stat;
2848 MPI_Recv(dfdx[0], natoms*DIM, mpi_type, node, node,
2849 cr->mpi_comm_mygroup, &stat);
2850 #undef mpi_type
2851 #endif
2854 row = (atom + node)*DIM + d;
2856 for (j = 0; j < natoms; j++)
2858 for (k = 0; k < DIM; k++)
2860 col = j*DIM + k;
2862 if (bSparse)
2864 if (col >= row && dfdx[j][k] != 0.0)
2866 gmx_sparsematrix_increment_value(sparse_matrix,
2867 row, col, dfdx[j][k]);
2870 else
2872 full_matrix[row*sz+col] = dfdx[j][k];
2879 if (bVerbose && fplog)
2881 fflush(fplog);
2884 /* write progress */
2885 if (MASTER(cr) && bVerbose)
2887 fprintf(stderr, "\rFinished step %d out of %d",
2888 min(atom+nnodes, natoms), natoms);
2889 fflush(stderr);
2893 if (MASTER(cr))
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);
2903 return 0;