Remove state duplication in serial runs
[gromacs.git] / src / gromacs / mdlib / minimize.cpp
blob1fee19c96766343ad6dddcd63f5c7b9aa27a7156
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,2015,2016, 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 /*! \internal \file
39 * \brief This file defines integrators for energy minimization
41 * \author Berk Hess <hess@kth.se>
42 * \author Erik Lindahl <erik@kth.se>
43 * \ingroup module_mdlib
45 #include "gmxpre.h"
47 #include "minimize.h"
49 #include "config.h"
51 #include <cmath>
52 #include <cstring>
53 #include <ctime>
55 #include <algorithm>
56 #include <vector>
58 #include "gromacs/commandline/filenm.h"
59 #include "gromacs/domdec/domdec.h"
60 #include "gromacs/domdec/domdec_struct.h"
61 #include "gromacs/ewald/pme.h"
62 #include "gromacs/fileio/confio.h"
63 #include "gromacs/fileio/mtxio.h"
64 #include "gromacs/gmxlib/network.h"
65 #include "gromacs/gmxlib/nrnb.h"
66 #include "gromacs/imd/imd.h"
67 #include "gromacs/linearalgebra/sparsematrix.h"
68 #include "gromacs/listed-forces/manage-threading.h"
69 #include "gromacs/math/functions.h"
70 #include "gromacs/math/vec.h"
71 #include "gromacs/mdlib/constr.h"
72 #include "gromacs/mdlib/force.h"
73 #include "gromacs/mdlib/forcerec.h"
74 #include "gromacs/mdlib/gmx_omp_nthreads.h"
75 #include "gromacs/mdlib/md_support.h"
76 #include "gromacs/mdlib/mdatoms.h"
77 #include "gromacs/mdlib/mdebin.h"
78 #include "gromacs/mdlib/mdrun.h"
79 #include "gromacs/mdlib/mdsetup.h"
80 #include "gromacs/mdlib/ns.h"
81 #include "gromacs/mdlib/shellfc.h"
82 #include "gromacs/mdlib/sim_util.h"
83 #include "gromacs/mdlib/tgroup.h"
84 #include "gromacs/mdlib/trajectory_writing.h"
85 #include "gromacs/mdlib/update.h"
86 #include "gromacs/mdlib/vsite.h"
87 #include "gromacs/mdtypes/commrec.h"
88 #include "gromacs/mdtypes/inputrec.h"
89 #include "gromacs/mdtypes/md_enums.h"
90 #include "gromacs/pbcutil/mshift.h"
91 #include "gromacs/pbcutil/pbc.h"
92 #include "gromacs/timing/wallcycle.h"
93 #include "gromacs/timing/walltime_accounting.h"
94 #include "gromacs/topology/mtop_util.h"
95 #include "gromacs/topology/topology.h"
96 #include "gromacs/utility/cstringutil.h"
97 #include "gromacs/utility/exceptions.h"
98 #include "gromacs/utility/fatalerror.h"
99 #include "gromacs/utility/logger.h"
100 #include "gromacs/utility/smalloc.h"
102 //! Utility structure for manipulating states during EM
103 typedef struct {
104 //! Copy of the global state
105 t_state s;
106 //! Force array
107 rvec *f;
108 //! Potential energy
109 real epot;
110 //! Norm of the force
111 real fnorm;
112 //! Maximum force
113 real fmax;
114 //! Direction
115 int a_fmax;
116 } em_state_t;
118 //! Initiate em_state_t structure and return pointer to it
119 static em_state_t *init_em_state()
121 em_state_t *ems;
123 snew(ems, 1);
125 /* does this need to be here? Should the array be declared differently (staticaly)in the state definition? */
126 snew(ems->s.lambda, efptNR);
128 return ems;
131 //! Print the EM starting conditions
132 static void print_em_start(FILE *fplog,
133 t_commrec *cr,
134 gmx_walltime_accounting_t walltime_accounting,
135 gmx_wallcycle_t wcycle,
136 const char *name)
138 walltime_accounting_start(walltime_accounting);
139 wallcycle_start(wcycle, ewcRUN);
140 print_start(fplog, cr, walltime_accounting, name);
143 //! Stop counting time for EM
144 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
145 gmx_wallcycle_t wcycle)
147 wallcycle_stop(wcycle, ewcRUN);
149 walltime_accounting_end(walltime_accounting);
152 //! Printing a log file and console header
153 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
155 fprintf(out, "\n");
156 fprintf(out, "%s:\n", minimizer);
157 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
158 fprintf(out, " Number of steps = %12d\n", nsteps);
161 //! Print warning message
162 static void warn_step(FILE *fp, real ftol, gmx_bool bLastStep, gmx_bool bConstrain)
164 char buffer[2048];
165 if (bLastStep)
167 sprintf(buffer,
168 "\nEnergy minimization reached the maximum number "
169 "of steps before the forces reached the requested "
170 "precision Fmax < %g.\n", ftol);
172 else
174 sprintf(buffer,
175 "\nEnergy minimization has stopped, but the forces have "
176 "not converged to the requested precision Fmax < %g (which "
177 "may not be possible for your system). It stopped "
178 "because the algorithm tried to make a new step whose size "
179 "was too small, or there was no change in the energy since "
180 "last step. Either way, we regard the minimization as "
181 "converged to within the available machine precision, "
182 "given your starting configuration and EM parameters.\n%s%s",
183 ftol,
184 sizeof(real) < sizeof(double) ?
185 "\nDouble precision normally gives you higher accuracy, but "
186 "this is often not needed for preparing to run molecular "
187 "dynamics.\n" :
189 bConstrain ?
190 "You might need to increase your constraint accuracy, or turn\n"
191 "off constraints altogether (set constraints = none in mdp file)\n" :
192 "");
194 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
197 //! Print message about convergence of the EM
198 static void print_converged(FILE *fp, const char *alg, real ftol,
199 gmx_int64_t count, gmx_bool bDone, gmx_int64_t nsteps,
200 real epot, real fmax, int nfmax, real fnorm)
202 char buf[STEPSTRSIZE];
204 if (bDone)
206 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
207 alg, ftol, gmx_step_str(count, buf));
209 else if (count < nsteps)
211 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
212 "but did not reach the requested Fmax < %g.\n",
213 alg, gmx_step_str(count, buf), ftol);
215 else
217 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
218 alg, ftol, gmx_step_str(count, buf));
221 #if GMX_DOUBLE
222 fprintf(fp, "Potential Energy = %21.14e\n", epot);
223 fprintf(fp, "Maximum force = %21.14e on atom %d\n", fmax, nfmax+1);
224 fprintf(fp, "Norm of force = %21.14e\n", fnorm);
225 #else
226 fprintf(fp, "Potential Energy = %14.7e\n", epot);
227 fprintf(fp, "Maximum force = %14.7e on atom %d\n", fmax, nfmax+1);
228 fprintf(fp, "Norm of force = %14.7e\n", fnorm);
229 #endif
232 //! Compute the norm and max of the force array in parallel
233 static void get_f_norm_max(t_commrec *cr,
234 t_grpopts *opts, t_mdatoms *mdatoms, rvec *f,
235 real *fnorm, real *fmax, int *a_fmax)
237 double fnorm2, *sum;
238 real fmax2, fam;
239 int la_max, a_max, start, end, i, m, gf;
241 /* This routine finds the largest force and returns it.
242 * On parallel machines the global max is taken.
244 fnorm2 = 0;
245 fmax2 = 0;
246 la_max = -1;
247 start = 0;
248 end = mdatoms->homenr;
249 if (mdatoms->cFREEZE)
251 for (i = start; i < end; i++)
253 gf = mdatoms->cFREEZE[i];
254 fam = 0;
255 for (m = 0; m < DIM; m++)
257 if (!opts->nFreeze[gf][m])
259 fam += gmx::square(f[i][m]);
262 fnorm2 += fam;
263 if (fam > fmax2)
265 fmax2 = fam;
266 la_max = i;
270 else
272 for (i = start; i < end; i++)
274 fam = norm2(f[i]);
275 fnorm2 += fam;
276 if (fam > fmax2)
278 fmax2 = fam;
279 la_max = i;
284 if (la_max >= 0 && DOMAINDECOMP(cr))
286 a_max = cr->dd->gatindex[la_max];
288 else
290 a_max = la_max;
292 if (PAR(cr))
294 snew(sum, 2*cr->nnodes+1);
295 sum[2*cr->nodeid] = fmax2;
296 sum[2*cr->nodeid+1] = a_max;
297 sum[2*cr->nnodes] = fnorm2;
298 gmx_sumd(2*cr->nnodes+1, sum, cr);
299 fnorm2 = sum[2*cr->nnodes];
300 /* Determine the global maximum */
301 for (i = 0; i < cr->nnodes; i++)
303 if (sum[2*i] > fmax2)
305 fmax2 = sum[2*i];
306 a_max = (int)(sum[2*i+1] + 0.5);
309 sfree(sum);
312 if (fnorm)
314 *fnorm = sqrt(fnorm2);
316 if (fmax)
318 *fmax = sqrt(fmax2);
320 if (a_fmax)
322 *a_fmax = a_max;
326 //! Compute the norm of the force
327 static void get_state_f_norm_max(t_commrec *cr,
328 t_grpopts *opts, t_mdatoms *mdatoms,
329 em_state_t *ems)
331 get_f_norm_max(cr, opts, mdatoms, ems->f, &ems->fnorm, &ems->fmax, &ems->a_fmax);
334 //! Initialize the energy minimization
335 void init_em(FILE *fplog, const char *title,
336 t_commrec *cr, t_inputrec *ir,
337 t_state *state_global, gmx_mtop_t *top_global,
338 em_state_t *ems, gmx_localtop_t **top,
339 rvec **f,
340 t_nrnb *nrnb, rvec mu_tot,
341 t_forcerec *fr, gmx_enerdata_t **enerd,
342 t_graph **graph, t_mdatoms *mdatoms, gmx_global_stat_t *gstat,
343 gmx_vsite_t *vsite, gmx_constr_t constr, gmx_shellfc_t **shellfc,
344 int nfile, const t_filenm fnm[],
345 gmx_mdoutf_t *outf, t_mdebin **mdebin,
346 int imdport, unsigned long gmx_unused Flags,
347 gmx_wallcycle_t wcycle)
349 int i;
350 real dvdl_constr;
352 if (fplog)
354 fprintf(fplog, "Initiating %s\n", title);
357 state_global->ngtc = 0;
359 /* Initialize lambda variables */
360 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, NULL);
362 init_nrnb(nrnb);
364 /* Interactive molecular dynamics */
365 init_IMD(ir, cr, top_global, fplog, 1, state_global->x,
366 nfile, fnm, NULL, imdport, Flags);
368 if (ir->eI == eiNM)
370 GMX_ASSERT(shellfc != NULL, "With NM we always support shells");
372 *shellfc = init_shell_flexcon(stdout,
373 top_global,
374 n_flexible_constraints(constr),
375 ir->nstcalcenergy,
376 DOMAINDECOMP(cr));
378 else
380 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI), "This else currently only handles energy minimizers, consider if your algorithm needs shell/flexible-constraint support");
382 /* With energy minimization, shells and flexible constraints are
383 * automatically minimized when treated like normal DOFS.
385 if (shellfc != NULL)
387 *shellfc = NULL;
391 if (DOMAINDECOMP(cr))
393 *top = dd_init_local_top(top_global);
395 dd_init_local_state(cr->dd, state_global, &ems->s);
397 *f = NULL;
399 /* Distribute the charge groups over the nodes from the master node */
400 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
401 state_global, top_global, ir,
402 &ems->s, &ems->f, mdatoms, *top,
403 fr, vsite, constr,
404 nrnb, NULL, FALSE);
405 dd_store_state(cr->dd, &ems->s);
407 *graph = NULL;
409 else
411 snew(*f, top_global->natoms);
413 /* Just copy the state */
414 ems->s = *state_global;
415 /* We need to allocate one element extra, since we might use
416 * (unaligned) 4-wide SIMD loads to access rvec entries.
418 snew(ems->s.x, ems->s.nalloc + 1);
419 snew(ems->f, ems->s.nalloc+1);
420 snew(ems->s.v, ems->s.nalloc+1);
421 for (i = 0; i < state_global->natoms; i++)
423 copy_rvec(state_global->x[i], ems->s.x[i]);
425 copy_mat(state_global->box, ems->s.box);
428 snew(*top, 1);
429 mdAlgorithmsSetupAtomData(cr, ir, top_global, *top, fr,
430 graph, mdatoms,
431 vsite, shellfc ? *shellfc : NULL);
433 update_mdatoms(mdatoms, state_global->lambda[efptFEP]);
435 if (vsite)
437 set_vsite_top(vsite, *top, mdatoms, cr);
441 if (constr)
443 if (ir->eConstrAlg == econtSHAKE &&
444 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
446 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
447 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
450 if (!DOMAINDECOMP(cr))
452 set_constraints(constr, *top, ir, mdatoms, cr);
455 if (!ir->bContinuation)
457 /* Constrain the starting coordinates */
458 dvdl_constr = 0;
459 constrain(PAR(cr) ? NULL : fplog, TRUE, TRUE, constr, &(*top)->idef,
460 ir, cr, -1, 0, 1.0, mdatoms,
461 ems->s.x, ems->s.x, NULL, fr->bMolPBC, ems->s.box,
462 ems->s.lambda[efptFEP], &dvdl_constr,
463 NULL, NULL, nrnb, econqCoord);
467 if (PAR(cr))
469 *gstat = global_stat_init(ir);
471 else
473 *gstat = NULL;
476 *outf = init_mdoutf(fplog, nfile, fnm, 0, cr, ir, top_global, NULL, wcycle);
478 snew(*enerd, 1);
479 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
480 *enerd);
482 if (mdebin != NULL)
484 /* Init bin for energy stuff */
485 *mdebin = init_mdebin(mdoutf_get_fp_ene(*outf), top_global, ir, NULL);
488 clear_rvec(mu_tot);
489 calc_shifts(ems->s.box, fr->shift_vec);
492 //! Finalize the minimization
493 static void finish_em(t_commrec *cr, gmx_mdoutf_t outf,
494 gmx_walltime_accounting_t walltime_accounting,
495 gmx_wallcycle_t wcycle)
497 if (!(cr->duty & DUTY_PME))
499 /* Tell the PME only node to finish */
500 gmx_pme_send_finish(cr);
503 done_mdoutf(outf);
505 em_time_end(walltime_accounting, wcycle);
508 //! Swap two different EM states during minimization
509 static void swap_em_state(em_state_t *ems1, em_state_t *ems2)
511 em_state_t tmp;
513 tmp = *ems1;
514 *ems1 = *ems2;
515 *ems2 = tmp;
518 //! Save the EM trajectory
519 static void write_em_traj(FILE *fplog, t_commrec *cr,
520 gmx_mdoutf_t outf,
521 gmx_bool bX, gmx_bool bF, const char *confout,
522 gmx_mtop_t *top_global,
523 t_inputrec *ir, gmx_int64_t step,
524 em_state_t *state,
525 t_state *state_global)
527 int mdof_flags = 0;
529 if (bX)
531 mdof_flags |= MDOF_X;
533 if (bF)
535 mdof_flags |= MDOF_F;
538 /* If we want IMD output, set appropriate MDOF flag */
539 if (ir->bIMD)
541 mdof_flags |= MDOF_IMD;
544 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
545 top_global, step, (double)step,
546 &state->s, state_global, state->f);
548 if (confout != NULL && MASTER(cr))
550 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
552 /* Make molecules whole only for confout writing */
553 do_pbc_mtop(fplog, ir->ePBC, state_global->box, top_global,
554 state_global->x);
557 write_sto_conf_mtop(confout,
558 *top_global->name, top_global,
559 state_global->x, NULL, ir->ePBC, state_global->box);
563 //! \brief Do one minimization step
565 // \returns true when the step succeeded, false when a constraint error occurred
566 static bool do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
567 gmx_bool bMolPBC,
568 em_state_t *ems1, real a, rvec *f, em_state_t *ems2,
569 gmx_constr_t constr, gmx_localtop_t *top,
570 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
571 gmx_int64_t count)
574 t_state *s1, *s2;
575 int start, end;
576 real dvdl_constr;
577 int nthreads gmx_unused;
579 bool validStep = true;
581 s1 = &ems1->s;
582 s2 = &ems2->s;
584 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
586 gmx_incons("state mismatch in do_em_step");
589 s2->flags = s1->flags;
591 if (s2->nalloc != s1->nalloc)
593 s2->nalloc = s1->nalloc;
594 /* We need to allocate one element extra, since we might use
595 * (unaligned) 4-wide SIMD loads to access rvec entries.
597 srenew(s2->x, s1->nalloc + 1);
598 srenew(ems2->f, s1->nalloc);
599 if (s2->flags & (1<<estCGP))
601 srenew(s2->cg_p, s1->nalloc + 1);
605 s2->natoms = s1->natoms;
606 copy_mat(s1->box, s2->box);
607 /* Copy free energy state */
608 for (int i = 0; i < efptNR; i++)
610 s2->lambda[i] = s1->lambda[i];
612 copy_mat(s1->box, s2->box);
614 start = 0;
615 end = md->homenr;
617 // cppcheck-suppress unreadVariable
618 nthreads = gmx_omp_nthreads_get(emntUpdate);
619 #pragma omp parallel num_threads(nthreads)
621 rvec *x1 = s1->x;
622 rvec *x2 = s2->x;
624 int gf = 0;
625 #pragma omp for schedule(static) nowait
626 for (int i = start; i < end; i++)
630 if (md->cFREEZE)
632 gf = md->cFREEZE[i];
634 for (int m = 0; m < DIM; m++)
636 if (ir->opts.nFreeze[gf][m])
638 x2[i][m] = x1[i][m];
640 else
642 x2[i][m] = x1[i][m] + a*f[i][m];
646 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
649 if (s2->flags & (1<<estCGP))
651 /* Copy the CG p vector */
652 rvec *p1 = s1->cg_p;
653 rvec *p2 = s2->cg_p;
654 #pragma omp for schedule(static) nowait
655 for (int i = start; i < end; i++)
657 // Trivial OpenMP block that does not throw
658 copy_rvec(p1[i], p2[i]);
662 if (DOMAINDECOMP(cr))
664 s2->ddp_count = s1->ddp_count;
665 if (s2->cg_gl_nalloc < s1->cg_gl_nalloc)
667 #pragma omp barrier
668 s2->cg_gl_nalloc = s1->cg_gl_nalloc;
671 /* We need to allocate one element extra, since we might use
672 * (unaligned) 4-wide SIMD loads to access rvec entries.
674 srenew(s2->cg_gl, s2->cg_gl_nalloc + 1);
676 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
677 #pragma omp barrier
679 s2->ncg_gl = s1->ncg_gl;
680 #pragma omp for schedule(static) nowait
681 for (int i = 0; i < s2->ncg_gl; i++)
683 s2->cg_gl[i] = s1->cg_gl[i];
685 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
689 if (constr)
691 wallcycle_start(wcycle, ewcCONSTR);
692 dvdl_constr = 0;
693 validStep =
694 constrain(NULL, TRUE, TRUE, constr, &top->idef,
695 ir, cr, count, 0, 1.0, md,
696 s1->x, s2->x, NULL, bMolPBC, s2->box,
697 s2->lambda[efptBONDED], &dvdl_constr,
698 NULL, NULL, nrnb, econqCoord);
699 wallcycle_stop(wcycle, ewcCONSTR);
701 // We should move this check to the different minimizers
702 if (!validStep && ir->eI != eiSteep)
704 gmx_fatal(FARGS, "The coordinates could not be constrained. Minimizer '%s' can not handle constraint failures, use minimizer '%s' before using '%s'.",
705 EI(ir->eI), EI(eiSteep), EI(ir->eI));
709 return validStep;
712 //! Prepare EM for using domain decomposition parallellization
713 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
714 gmx_mtop_t *top_global, t_inputrec *ir,
715 em_state_t *ems, gmx_localtop_t *top,
716 t_mdatoms *mdatoms, t_forcerec *fr,
717 gmx_vsite_t *vsite, gmx_constr_t constr,
718 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
720 /* Repartition the domain decomposition */
721 dd_partition_system(fplog, step, cr, FALSE, 1,
722 NULL, top_global, ir,
723 &ems->s, &ems->f,
724 mdatoms, top, fr, vsite, constr,
725 nrnb, wcycle, FALSE);
726 dd_store_state(cr->dd, &ems->s);
729 //! De one energy evaluation
730 static void evaluate_energy(FILE *fplog, t_commrec *cr,
731 gmx_mtop_t *top_global,
732 em_state_t *ems, gmx_localtop_t *top,
733 t_inputrec *inputrec,
734 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
735 gmx_global_stat_t gstat,
736 gmx_vsite_t *vsite, gmx_constr_t constr,
737 t_fcdata *fcd,
738 t_graph *graph, t_mdatoms *mdatoms,
739 t_forcerec *fr, rvec mu_tot,
740 gmx_enerdata_t *enerd, tensor vir, tensor pres,
741 gmx_int64_t count, gmx_bool bFirst)
743 real t;
744 gmx_bool bNS;
745 tensor force_vir, shake_vir, ekin;
746 real dvdl_constr, prescorr, enercorr, dvdlcorr;
747 real terminate = 0;
749 /* Set the time to the initial time, the time does not change during EM */
750 t = inputrec->init_t;
752 if (bFirst ||
753 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
755 /* This is the first state or an old state used before the last ns */
756 bNS = TRUE;
758 else
760 bNS = FALSE;
761 if (inputrec->nstlist > 0)
763 bNS = TRUE;
767 if (vsite)
769 construct_vsites(vsite, ems->s.x, 1, NULL,
770 top->idef.iparams, top->idef.il,
771 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
774 if (DOMAINDECOMP(cr) && bNS)
776 /* Repartition the domain decomposition */
777 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
778 ems, top, mdatoms, fr, vsite, constr,
779 nrnb, wcycle);
782 /* Calc force & energy on new trial position */
783 /* do_force always puts the charge groups in the box and shifts again
784 * We do not unshift, so molecules are always whole in congrad.c
786 do_force(fplog, cr, inputrec,
787 count, nrnb, wcycle, top, &top_global->groups,
788 ems->s.box, ems->s.x, &ems->s.hist,
789 ems->f, force_vir, mdatoms, enerd, fcd,
790 ems->s.lambda, graph, fr, vsite, mu_tot, t, NULL, NULL, TRUE,
791 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
792 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
793 (bNS ? GMX_FORCE_NS : 0));
795 /* Clear the unused shake virial and pressure */
796 clear_mat(shake_vir);
797 clear_mat(pres);
799 /* Communicate stuff when parallel */
800 if (PAR(cr) && inputrec->eI != eiNM)
802 wallcycle_start(wcycle, ewcMoveE);
804 global_stat(gstat, cr, enerd, force_vir, shake_vir, mu_tot,
805 inputrec, NULL, NULL, NULL, 1, &terminate,
806 NULL, FALSE,
807 CGLO_ENERGY |
808 CGLO_PRESSURE |
809 CGLO_CONSTRAINT);
811 wallcycle_stop(wcycle, ewcMoveE);
814 /* Calculate long range corrections to pressure and energy */
815 calc_dispcorr(inputrec, fr, ems->s.box, ems->s.lambda[efptVDW],
816 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
817 enerd->term[F_DISPCORR] = enercorr;
818 enerd->term[F_EPOT] += enercorr;
819 enerd->term[F_PRES] += prescorr;
820 enerd->term[F_DVDL] += dvdlcorr;
822 ems->epot = enerd->term[F_EPOT];
824 if (constr)
826 /* Project out the constraint components of the force */
827 wallcycle_start(wcycle, ewcCONSTR);
828 dvdl_constr = 0;
829 constrain(NULL, FALSE, FALSE, constr, &top->idef,
830 inputrec, cr, count, 0, 1.0, mdatoms,
831 ems->s.x, ems->f, ems->f, fr->bMolPBC, ems->s.box,
832 ems->s.lambda[efptBONDED], &dvdl_constr,
833 NULL, &shake_vir, nrnb, econqForceDispl);
834 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
835 m_add(force_vir, shake_vir, vir);
836 wallcycle_stop(wcycle, ewcCONSTR);
838 else
840 copy_mat(force_vir, vir);
843 clear_mat(ekin);
844 enerd->term[F_PRES] =
845 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
847 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
849 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
851 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
855 //! Parallel utility summing energies and forces
856 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
857 gmx_mtop_t *top_global,
858 em_state_t *s_min, em_state_t *s_b)
860 rvec *fm, *fb, *fmg;
861 t_block *cgs_gl;
862 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
863 double partsum;
864 unsigned char *grpnrFREEZE;
866 if (debug)
868 fprintf(debug, "Doing reorder_partsum\n");
871 fm = s_min->f;
872 fb = s_b->f;
874 cgs_gl = dd_charge_groups_global(cr->dd);
875 index = cgs_gl->index;
877 /* Collect fm in a global vector fmg.
878 * This conflicts with the spirit of domain decomposition,
879 * but to fully optimize this a much more complicated algorithm is required.
881 snew(fmg, top_global->natoms);
883 ncg = s_min->s.ncg_gl;
884 cg_gl = s_min->s.cg_gl;
885 i = 0;
886 for (c = 0; c < ncg; c++)
888 cg = cg_gl[c];
889 a0 = index[cg];
890 a1 = index[cg+1];
891 for (a = a0; a < a1; a++)
893 copy_rvec(fm[i], fmg[a]);
894 i++;
897 gmx_sum(top_global->natoms*3, fmg[0], cr);
899 /* Now we will determine the part of the sum for the cgs in state s_b */
900 ncg = s_b->s.ncg_gl;
901 cg_gl = s_b->s.cg_gl;
902 partsum = 0;
903 i = 0;
904 gf = 0;
905 grpnrFREEZE = top_global->groups.grpnr[egcFREEZE];
906 for (c = 0; c < ncg; c++)
908 cg = cg_gl[c];
909 a0 = index[cg];
910 a1 = index[cg+1];
911 for (a = a0; a < a1; a++)
913 if (mdatoms->cFREEZE && grpnrFREEZE)
915 gf = grpnrFREEZE[i];
917 for (m = 0; m < DIM; m++)
919 if (!opts->nFreeze[gf][m])
921 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
924 i++;
928 sfree(fmg);
930 return partsum;
933 //! Print some stuff, like beta, whatever that means.
934 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
935 gmx_mtop_t *top_global,
936 em_state_t *s_min, em_state_t *s_b)
938 rvec *fm, *fb;
939 double sum;
940 int gf, i, m;
942 /* This is just the classical Polak-Ribiere calculation of beta;
943 * it looks a bit complicated since we take freeze groups into account,
944 * and might have to sum it in parallel runs.
947 if (!DOMAINDECOMP(cr) ||
948 (s_min->s.ddp_count == cr->dd->ddp_count &&
949 s_b->s.ddp_count == cr->dd->ddp_count))
951 fm = s_min->f;
952 fb = s_b->f;
953 sum = 0;
954 gf = 0;
955 /* This part of code can be incorrect with DD,
956 * since the atom ordering in s_b and s_min might differ.
958 for (i = 0; i < mdatoms->homenr; i++)
960 if (mdatoms->cFREEZE)
962 gf = mdatoms->cFREEZE[i];
964 for (m = 0; m < DIM; m++)
966 if (!opts->nFreeze[gf][m])
968 sum += (fb[i][m] - fm[i][m])*fb[i][m];
973 else
975 /* We need to reorder cgs while summing */
976 sum = reorder_partsum(cr, opts, mdatoms, top_global, s_min, s_b);
978 if (PAR(cr))
980 gmx_sumd(1, &sum, cr);
983 return sum/gmx::square(s_min->fnorm);
986 namespace gmx
989 /*! \brief Do conjugate gradients minimization
990 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
991 int nfile, const t_filenm fnm[],
992 const gmx_output_env_t *oenv, gmx_bool bVerbose,
993 int nstglobalcomm,
994 gmx_vsite_t *vsite, gmx_constr_t constr,
995 int stepout,
996 t_inputrec *inputrec,
997 gmx_mtop_t *top_global, t_fcdata *fcd,
998 t_state *state_global,
999 t_mdatoms *mdatoms,
1000 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1001 gmx_edsam_t ed,
1002 t_forcerec *fr,
1003 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
1004 gmx_membed_t gmx_unused *membed,
1005 real cpt_period, real max_hours,
1006 int imdport,
1007 unsigned long Flags,
1008 gmx_walltime_accounting_t walltime_accounting)
1010 double do_cg(FILE *fplog, t_commrec *cr, const gmx::MDLogger gmx_unused &mdlog,
1011 int nfile, const t_filenm fnm[],
1012 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
1013 int gmx_unused nstglobalcomm,
1014 gmx_vsite_t *vsite, gmx_constr_t constr,
1015 int gmx_unused stepout,
1016 t_inputrec *inputrec,
1017 gmx_mtop_t *top_global, t_fcdata *fcd,
1018 t_state *state_global,
1019 t_mdatoms *mdatoms,
1020 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1021 gmx_edsam_t gmx_unused ed,
1022 t_forcerec *fr,
1023 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1024 gmx_membed_t gmx_unused *membed,
1025 real gmx_unused cpt_period, real gmx_unused max_hours,
1026 int imdport,
1027 unsigned long gmx_unused Flags,
1028 gmx_walltime_accounting_t walltime_accounting)
1030 const char *CG = "Polak-Ribiere Conjugate Gradients";
1032 em_state_t *s_min, *s_a, *s_b, *s_c;
1033 gmx_localtop_t *top;
1034 gmx_enerdata_t *enerd;
1035 rvec *f;
1036 gmx_global_stat_t gstat;
1037 t_graph *graph;
1038 rvec *p, *sf;
1039 double gpa, gpb, gpc, tmp, minstep;
1040 real fnormn;
1041 real stepsize;
1042 real a, b, c, beta = 0.0;
1043 real epot_repl = 0;
1044 real pnorm;
1045 t_mdebin *mdebin;
1046 gmx_bool converged, foundlower;
1047 rvec mu_tot;
1048 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1049 tensor vir, pres;
1050 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1051 gmx_mdoutf_t outf;
1052 int i, m, gf, step, nminstep;
1054 step = 0;
1056 s_min = init_em_state();
1057 s_a = init_em_state();
1058 s_b = init_em_state();
1059 s_c = init_em_state();
1061 /* Init em and store the local state in s_min */
1062 init_em(fplog, CG, cr, inputrec,
1063 state_global, top_global, s_min, &top, &f,
1064 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat,
1065 vsite, constr, NULL,
1066 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
1068 /* Print to log file */
1069 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1071 /* Max number of steps */
1072 number_steps = inputrec->nsteps;
1074 if (MASTER(cr))
1076 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1078 if (fplog)
1080 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1083 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1084 /* do_force always puts the charge groups in the box and shifts again
1085 * We do not unshift, so molecules are always whole in congrad.c
1087 evaluate_energy(fplog, cr,
1088 top_global, s_min, top,
1089 inputrec, nrnb, wcycle, gstat,
1090 vsite, constr, fcd, graph, mdatoms, fr,
1091 mu_tot, enerd, vir, pres, -1, TRUE);
1092 where();
1094 if (MASTER(cr))
1096 /* Copy stuff to the energy bin for easy printing etc. */
1097 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1098 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1099 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1101 print_ebin_header(fplog, step, step);
1102 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1103 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1105 where();
1107 /* Estimate/guess the initial stepsize */
1108 stepsize = inputrec->em_stepsize/s_min->fnorm;
1110 if (MASTER(cr))
1112 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1113 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1114 s_min->fmax, s_min->a_fmax+1);
1115 fprintf(stderr, " F-Norm = %12.5e\n",
1116 s_min->fnorm/sqrtNumAtoms);
1117 fprintf(stderr, "\n");
1118 /* and copy to the log file too... */
1119 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1120 s_min->fmax, s_min->a_fmax+1);
1121 fprintf(fplog, " F-Norm = %12.5e\n",
1122 s_min->fnorm/sqrtNumAtoms);
1123 fprintf(fplog, "\n");
1125 /* Start the loop over CG steps.
1126 * Each successful step is counted, and we continue until
1127 * we either converge or reach the max number of steps.
1129 converged = FALSE;
1130 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1133 /* start taking steps in a new direction
1134 * First time we enter the routine, beta=0, and the direction is
1135 * simply the negative gradient.
1138 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1139 p = s_min->s.cg_p;
1140 sf = s_min->f;
1141 gpa = 0;
1142 gf = 0;
1143 for (i = 0; i < mdatoms->homenr; i++)
1145 if (mdatoms->cFREEZE)
1147 gf = mdatoms->cFREEZE[i];
1149 for (m = 0; m < DIM; m++)
1151 if (!inputrec->opts.nFreeze[gf][m])
1153 p[i][m] = sf[i][m] + beta*p[i][m];
1154 gpa -= p[i][m]*sf[i][m];
1155 /* f is negative gradient, thus the sign */
1157 else
1159 p[i][m] = 0;
1164 /* Sum the gradient along the line across CPUs */
1165 if (PAR(cr))
1167 gmx_sumd(1, &gpa, cr);
1170 /* Calculate the norm of the search vector */
1171 get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1173 /* Just in case stepsize reaches zero due to numerical precision... */
1174 if (stepsize <= 0)
1176 stepsize = inputrec->em_stepsize/pnorm;
1180 * Double check the value of the derivative in the search direction.
1181 * If it is positive it must be due to the old information in the
1182 * CG formula, so just remove that and start over with beta=0.
1183 * This corresponds to a steepest descent step.
1185 if (gpa > 0)
1187 beta = 0;
1188 step--; /* Don't count this step since we are restarting */
1189 continue; /* Go back to the beginning of the big for-loop */
1192 /* Calculate minimum allowed stepsize, before the average (norm)
1193 * relative change in coordinate is smaller than precision
1195 minstep = 0;
1196 for (i = 0; i < mdatoms->homenr; i++)
1198 for (m = 0; m < DIM; m++)
1200 tmp = fabs(s_min->s.x[i][m]);
1201 if (tmp < 1.0)
1203 tmp = 1.0;
1205 tmp = p[i][m]/tmp;
1206 minstep += tmp*tmp;
1209 /* Add up from all CPUs */
1210 if (PAR(cr))
1212 gmx_sumd(1, &minstep, cr);
1215 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1217 if (stepsize < minstep)
1219 converged = TRUE;
1220 break;
1223 /* Write coordinates if necessary */
1224 do_x = do_per_step(step, inputrec->nstxout);
1225 do_f = do_per_step(step, inputrec->nstfout);
1227 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1228 top_global, inputrec, step,
1229 s_min, state_global);
1231 /* Take a step downhill.
1232 * In theory, we should minimize the function along this direction.
1233 * That is quite possible, but it turns out to take 5-10 function evaluations
1234 * for each line. However, we dont really need to find the exact minimum -
1235 * it is much better to start a new CG step in a modified direction as soon
1236 * as we are close to it. This will save a lot of energy evaluations.
1238 * In practice, we just try to take a single step.
1239 * If it worked (i.e. lowered the energy), we increase the stepsize but
1240 * the continue straight to the next CG step without trying to find any minimum.
1241 * If it didn't work (higher energy), there must be a minimum somewhere between
1242 * the old position and the new one.
1244 * Due to the finite numerical accuracy, it turns out that it is a good idea
1245 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1246 * This leads to lower final energies in the tests I've done. / Erik
1248 s_a->epot = s_min->epot;
1249 a = 0.0;
1250 c = a + stepsize; /* reference position along line is zero */
1252 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1254 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1255 s_min, top, mdatoms, fr, vsite, constr,
1256 nrnb, wcycle);
1259 /* Take a trial step (new coords in s_c) */
1260 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1261 constr, top, nrnb, wcycle, -1);
1263 neval++;
1264 /* Calculate energy for the trial step */
1265 evaluate_energy(fplog, cr,
1266 top_global, s_c, top,
1267 inputrec, nrnb, wcycle, gstat,
1268 vsite, constr, fcd, graph, mdatoms, fr,
1269 mu_tot, enerd, vir, pres, -1, FALSE);
1271 /* Calc derivative along line */
1272 p = s_c->s.cg_p;
1273 sf = s_c->f;
1274 gpc = 0;
1275 for (i = 0; i < mdatoms->homenr; i++)
1277 for (m = 0; m < DIM; m++)
1279 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1282 /* Sum the gradient along the line across CPUs */
1283 if (PAR(cr))
1285 gmx_sumd(1, &gpc, cr);
1288 /* This is the max amount of increase in energy we tolerate */
1289 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1291 /* Accept the step if the energy is lower, or if it is not significantly higher
1292 * and the line derivative is still negative.
1294 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1296 foundlower = TRUE;
1297 /* Great, we found a better energy. Increase step for next iteration
1298 * if we are still going down, decrease it otherwise
1300 if (gpc < 0)
1302 stepsize *= 1.618034; /* The golden section */
1304 else
1306 stepsize *= 0.618034; /* 1/golden section */
1309 else
1311 /* New energy is the same or higher. We will have to do some work
1312 * to find a smaller value in the interval. Take smaller step next time!
1314 foundlower = FALSE;
1315 stepsize *= 0.618034;
1321 /* OK, if we didn't find a lower value we will have to locate one now - there must
1322 * be one in the interval [a=0,c].
1323 * The same thing is valid here, though: Don't spend dozens of iterations to find
1324 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1325 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1327 * I also have a safeguard for potentially really pathological functions so we never
1328 * take more than 20 steps before we give up ...
1330 * If we already found a lower value we just skip this step and continue to the update.
1332 if (!foundlower)
1334 nminstep = 0;
1338 /* Select a new trial point.
1339 * If the derivatives at points a & c have different sign we interpolate to zero,
1340 * otherwise just do a bisection.
1342 if (gpa < 0 && gpc > 0)
1344 b = a + gpa*(a-c)/(gpc-gpa);
1346 else
1348 b = 0.5*(a+c);
1351 /* safeguard if interpolation close to machine accuracy causes errors:
1352 * never go outside the interval
1354 if (b <= a || b >= c)
1356 b = 0.5*(a+c);
1359 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1361 /* Reload the old state */
1362 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1363 s_min, top, mdatoms, fr, vsite, constr,
1364 nrnb, wcycle);
1367 /* Take a trial step to this new point - new coords in s_b */
1368 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1369 constr, top, nrnb, wcycle, -1);
1371 neval++;
1372 /* Calculate energy for the trial step */
1373 evaluate_energy(fplog, cr,
1374 top_global, s_b, top,
1375 inputrec, nrnb, wcycle, gstat,
1376 vsite, constr, fcd, graph, mdatoms, fr,
1377 mu_tot, enerd, vir, pres, -1, FALSE);
1379 /* p does not change within a step, but since the domain decomposition
1380 * might change, we have to use cg_p of s_b here.
1382 p = s_b->s.cg_p;
1383 sf = s_b->f;
1384 gpb = 0;
1385 for (i = 0; i < mdatoms->homenr; i++)
1387 for (m = 0; m < DIM; m++)
1389 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1392 /* Sum the gradient along the line across CPUs */
1393 if (PAR(cr))
1395 gmx_sumd(1, &gpb, cr);
1398 if (debug)
1400 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1401 s_a->epot, s_b->epot, s_c->epot, gpb);
1404 epot_repl = s_b->epot;
1406 /* Keep one of the intervals based on the value of the derivative at the new point */
1407 if (gpb > 0)
1409 /* Replace c endpoint with b */
1410 swap_em_state(s_b, s_c);
1411 c = b;
1412 gpc = gpb;
1414 else
1416 /* Replace a endpoint with b */
1417 swap_em_state(s_b, s_a);
1418 a = b;
1419 gpa = gpb;
1423 * Stop search as soon as we find a value smaller than the endpoints.
1424 * Never run more than 20 steps, no matter what.
1426 nminstep++;
1428 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1429 (nminstep < 20));
1431 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1432 nminstep >= 20)
1434 /* OK. We couldn't find a significantly lower energy.
1435 * If beta==0 this was steepest descent, and then we give up.
1436 * If not, set beta=0 and restart with steepest descent before quitting.
1438 if (beta == 0.0)
1440 /* Converged */
1441 converged = TRUE;
1442 break;
1444 else
1446 /* Reset memory before giving up */
1447 beta = 0.0;
1448 continue;
1452 /* Select min energy state of A & C, put the best in B.
1454 if (s_c->epot < s_a->epot)
1456 if (debug)
1458 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1459 s_c->epot, s_a->epot);
1461 swap_em_state(s_b, s_c);
1462 gpb = gpc;
1464 else
1466 if (debug)
1468 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1469 s_a->epot, s_c->epot);
1471 swap_em_state(s_b, s_a);
1472 gpb = gpa;
1476 else
1478 if (debug)
1480 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1481 s_c->epot);
1483 swap_em_state(s_b, s_c);
1484 gpb = gpc;
1487 /* new search direction */
1488 /* beta = 0 means forget all memory and restart with steepest descents. */
1489 if (nstcg && ((step % nstcg) == 0))
1491 beta = 0.0;
1493 else
1495 /* s_min->fnorm cannot be zero, because then we would have converged
1496 * and broken out.
1499 /* Polak-Ribiere update.
1500 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1502 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1504 /* Limit beta to prevent oscillations */
1505 if (fabs(beta) > 5.0)
1507 beta = 0.0;
1511 /* update positions */
1512 swap_em_state(s_min, s_b);
1513 gpa = gpb;
1515 /* Print it if necessary */
1516 if (MASTER(cr))
1518 if (bVerbose)
1520 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1521 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1522 step, s_min->epot, s_min->fnorm/sqrtNumAtoms,
1523 s_min->fmax, s_min->a_fmax+1);
1524 fflush(stderr);
1526 /* Store the new (lower) energies */
1527 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1528 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1529 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1531 do_log = do_per_step(step, inputrec->nstlog);
1532 do_ene = do_per_step(step, inputrec->nstenergy);
1534 /* Prepare IMD energy record, if bIMD is TRUE. */
1535 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1537 if (do_log)
1539 print_ebin_header(fplog, step, step);
1541 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1542 do_log ? fplog : NULL, step, step, eprNORMAL,
1543 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1546 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1547 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
1549 IMD_send_positions(inputrec->imd);
1552 /* Stop when the maximum force lies below tolerance.
1553 * If we have reached machine precision, converged is already set to true.
1555 converged = converged || (s_min->fmax < inputrec->em_tol);
1557 } /* End of the loop */
1559 /* IMD cleanup, if bIMD is TRUE. */
1560 IMD_finalize(inputrec->bIMD, inputrec->imd);
1562 if (converged)
1564 step--; /* we never took that last step in this case */
1567 if (s_min->fmax > inputrec->em_tol)
1569 if (MASTER(cr))
1571 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1572 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1574 converged = FALSE;
1577 if (MASTER(cr))
1579 /* If we printed energy and/or logfile last step (which was the last step)
1580 * we don't have to do it again, but otherwise print the final values.
1582 if (!do_log)
1584 /* Write final value to log since we didn't do anything the last step */
1585 print_ebin_header(fplog, step, step);
1587 if (!do_ene || !do_log)
1589 /* Write final energy file entries */
1590 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1591 !do_log ? fplog : NULL, step, step, eprNORMAL,
1592 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1596 /* Print some stuff... */
1597 if (MASTER(cr))
1599 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1602 /* IMPORTANT!
1603 * For accurate normal mode calculation it is imperative that we
1604 * store the last conformation into the full precision binary trajectory.
1606 * However, we should only do it if we did NOT already write this step
1607 * above (which we did if do_x or do_f was true).
1609 do_x = !do_per_step(step, inputrec->nstxout);
1610 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1612 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1613 top_global, inputrec, step,
1614 s_min, state_global);
1617 if (MASTER(cr))
1619 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1620 fnormn = s_min->fnorm/sqrtNumAtoms;
1621 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1622 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1623 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1624 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1626 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1629 finish_em(cr, outf, walltime_accounting, wcycle);
1631 /* To print the actual number of steps we needed somewhere */
1632 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1634 return 0;
1635 } /* That's all folks */
1638 /*! \brief Do L-BFGS conjugate gradients minimization
1639 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
1640 int nfile, const t_filenm fnm[],
1641 const gmx_output_env_t *oenv, gmx_bool bVerbose,
1642 int nstglobalcomm,
1643 gmx_vsite_t *vsite, gmx_constr_t constr,
1644 int stepout,
1645 t_inputrec *inputrec,
1646 gmx_mtop_t *top_global, t_fcdata *fcd,
1647 t_state *state_global,
1648 t_mdatoms *mdatoms,
1649 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1650 gmx_edsam_t ed,
1651 t_forcerec *fr,
1652 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
1653 real cpt_period, real max_hours,
1654 int imdport,
1655 unsigned long Flags,
1656 gmx_walltime_accounting_t walltime_accounting)
1658 double do_lbfgs(FILE *fplog, t_commrec *cr, const gmx::MDLogger gmx_unused &mdlog,
1659 int nfile, const t_filenm fnm[],
1660 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
1661 int gmx_unused nstglobalcomm,
1662 gmx_vsite_t *vsite, gmx_constr_t constr,
1663 int gmx_unused stepout,
1664 t_inputrec *inputrec,
1665 gmx_mtop_t *top_global, t_fcdata *fcd,
1666 t_state *state_global,
1667 t_mdatoms *mdatoms,
1668 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1669 gmx_edsam_t gmx_unused ed,
1670 t_forcerec *fr,
1671 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1672 gmx_membed_t gmx_unused *membed,
1673 real gmx_unused cpt_period, real gmx_unused max_hours,
1674 int imdport,
1675 unsigned long gmx_unused Flags,
1676 gmx_walltime_accounting_t walltime_accounting)
1678 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1679 em_state_t ems;
1680 gmx_localtop_t *top;
1681 gmx_enerdata_t *enerd;
1682 rvec *f;
1683 gmx_global_stat_t gstat;
1684 t_graph *graph;
1685 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1686 double stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1687 real *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1688 real *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1689 real a, b, c, maxdelta, delta;
1690 real diag, Epot0, Epot, EpotA, EpotB, EpotC;
1691 real dgdx, dgdg, sq, yr, beta;
1692 t_mdebin *mdebin;
1693 gmx_bool converged;
1694 rvec mu_tot;
1695 real fnorm, fmax;
1696 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1697 tensor vir, pres;
1698 int start, end, number_steps;
1699 gmx_mdoutf_t outf;
1700 int i, k, m, n, nfmax, gf, step;
1701 int mdof_flags;
1703 if (PAR(cr))
1705 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1708 if (NULL != constr)
1710 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).");
1713 n = 3*state_global->natoms;
1714 nmaxcorr = inputrec->nbfgscorr;
1716 /* Allocate memory */
1717 /* Use pointers to real so we dont have to loop over both atoms and
1718 * dimensions all the time...
1719 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1720 * that point to the same memory.
1722 snew(xa, n);
1723 snew(xb, n);
1724 snew(xc, n);
1725 snew(fa, n);
1726 snew(fb, n);
1727 snew(fc, n);
1728 snew(frozen, n);
1730 snew(p, n);
1731 snew(lastx, n);
1732 snew(lastf, n);
1733 snew(rho, nmaxcorr);
1734 snew(alpha, nmaxcorr);
1736 snew(dx, nmaxcorr);
1737 for (i = 0; i < nmaxcorr; i++)
1739 snew(dx[i], n);
1742 snew(dg, nmaxcorr);
1743 for (i = 0; i < nmaxcorr; i++)
1745 snew(dg[i], n);
1748 step = 0;
1749 neval = 0;
1751 /* Init em */
1752 init_em(fplog, LBFGS, cr, inputrec,
1753 state_global, top_global, &ems, &top, &f,
1754 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat,
1755 vsite, constr, NULL,
1756 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
1757 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1758 * so we free some memory again.
1760 sfree(ems.s.x);
1761 sfree(ems.f);
1763 xx = (real *)state_global->x;
1764 ff = (real *)f;
1766 start = 0;
1767 end = mdatoms->homenr;
1769 /* Print to log file */
1770 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1772 do_log = do_ene = do_x = do_f = TRUE;
1774 /* Max number of steps */
1775 number_steps = inputrec->nsteps;
1777 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1778 gf = 0;
1779 for (i = start; i < end; i++)
1781 if (mdatoms->cFREEZE)
1783 gf = mdatoms->cFREEZE[i];
1785 for (m = 0; m < DIM; m++)
1787 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1790 if (MASTER(cr))
1792 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1794 if (fplog)
1796 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1799 if (vsite)
1801 construct_vsites(vsite, state_global->x, 1, NULL,
1802 top->idef.iparams, top->idef.il,
1803 fr->ePBC, fr->bMolPBC, cr, state_global->box);
1806 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1807 /* do_force always puts the charge groups in the box and shifts again
1808 * We do not unshift, so molecules are always whole
1810 neval++;
1811 ems.s.x = state_global->x;
1812 ems.f = f;
1813 evaluate_energy(fplog, cr,
1814 top_global, &ems, top,
1815 inputrec, nrnb, wcycle, gstat,
1816 vsite, constr, fcd, graph, mdatoms, fr,
1817 mu_tot, enerd, vir, pres, -1, TRUE);
1818 where();
1820 if (MASTER(cr))
1822 /* Copy stuff to the energy bin for easy printing etc. */
1823 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1824 mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
1825 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1827 print_ebin_header(fplog, step, step);
1828 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1829 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1831 where();
1833 /* This is the starting energy */
1834 Epot = enerd->term[F_EPOT];
1836 fnorm = ems.fnorm;
1837 fmax = ems.fmax;
1838 nfmax = ems.a_fmax;
1840 /* Set the initial step.
1841 * since it will be multiplied by the non-normalized search direction
1842 * vector (force vector the first time), we scale it by the
1843 * norm of the force.
1846 if (MASTER(cr))
1848 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1849 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1850 fprintf(stderr, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1851 fprintf(stderr, " F-Norm = %12.5e\n", fnorm/sqrtNumAtoms);
1852 fprintf(stderr, "\n");
1853 /* and copy to the log file too... */
1854 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1855 fprintf(fplog, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1856 fprintf(fplog, " F-Norm = %12.5e\n", fnorm/sqrtNumAtoms);
1857 fprintf(fplog, "\n");
1860 // Point is an index to the memory of search directions, where 0 is the first one.
1861 point = 0;
1863 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1864 for (i = 0; i < n; i++)
1866 if (!frozen[i])
1868 dx[point][i] = ff[i]; /* Initial search direction */
1870 else
1872 dx[point][i] = 0;
1876 // Stepsize will be modified during the search, and actually it is not critical
1877 // (the main efficiency in the algorithm comes from changing directions), but
1878 // we still need an initial value, so estimate it as the inverse of the norm
1879 // so we take small steps where the potential fluctuates a lot.
1880 stepsize = 1.0/fnorm;
1882 /* Start the loop over BFGS steps.
1883 * Each successful step is counted, and we continue until
1884 * we either converge or reach the max number of steps.
1887 ncorr = 0;
1889 /* Set the gradient from the force */
1890 converged = FALSE;
1891 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1894 /* Write coordinates if necessary */
1895 do_x = do_per_step(step, inputrec->nstxout);
1896 do_f = do_per_step(step, inputrec->nstfout);
1898 mdof_flags = 0;
1899 if (do_x)
1901 mdof_flags |= MDOF_X;
1904 if (do_f)
1906 mdof_flags |= MDOF_F;
1909 if (inputrec->bIMD)
1911 mdof_flags |= MDOF_IMD;
1914 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1915 top_global, step, (real)step, state_global, state_global, f);
1917 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1919 /* make s a pointer to current search direction - point=0 first time we get here */
1920 s = dx[point];
1922 // calculate line gradient in position A
1923 for (gpa = 0, i = 0; i < n; i++)
1925 gpa -= s[i]*ff[i];
1928 /* Calculate minimum allowed stepsize along the line, before the average (norm)
1929 * relative change in coordinate is smaller than precision
1931 for (minstep = 0, i = 0; i < n; i++)
1933 tmp = fabs(xx[i]);
1934 if (tmp < 1.0)
1936 tmp = 1.0;
1938 tmp = s[i]/tmp;
1939 minstep += tmp*tmp;
1941 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1943 if (stepsize < minstep)
1945 converged = TRUE;
1946 break;
1949 // Before taking any steps along the line, store the old position
1950 for (i = 0; i < n; i++)
1952 lastx[i] = xx[i];
1953 lastf[i] = ff[i];
1955 Epot0 = Epot;
1957 for (i = 0; i < n; i++)
1959 xa[i] = xx[i];
1962 /* Take a step downhill.
1963 * In theory, we should find the actual minimum of the function in this
1964 * direction, somewhere along the line.
1965 * That is quite possible, but it turns out to take 5-10 function evaluations
1966 * for each line. However, we dont really need to find the exact minimum -
1967 * it is much better to start a new BFGS step in a modified direction as soon
1968 * as we are close to it. This will save a lot of energy evaluations.
1970 * In practice, we just try to take a single step.
1971 * If it worked (i.e. lowered the energy), we increase the stepsize but
1972 * continue straight to the next BFGS step without trying to find any minimum,
1973 * i.e. we change the search direction too. If the line was smooth, it is
1974 * likely we are in a smooth region, and then it makes sense to take longer
1975 * steps in the modified search direction too.
1977 * If it didn't work (higher energy), there must be a minimum somewhere between
1978 * the old position and the new one. Then we need to start by finding a lower
1979 * value before we change search direction. Since the energy was apparently
1980 * quite rough, we need to decrease the step size.
1982 * Due to the finite numerical accuracy, it turns out that it is a good idea
1983 * to accept a SMALL increase in energy, if the derivative is still downhill.
1984 * This leads to lower final energies in the tests I've done. / Erik
1987 // State "A" is the first position along the line.
1988 // reference position along line is initially zero
1989 EpotA = Epot0;
1990 a = 0.0;
1992 // Check stepsize first. We do not allow displacements
1993 // larger than emstep.
1997 // Pick a new position C by adding stepsize to A.
1998 c = a + stepsize;
2000 // Calculate what the largest change in any individual coordinate
2001 // would be (translation along line * gradient along line)
2002 maxdelta = 0;
2003 for (i = 0; i < n; i++)
2005 delta = c*s[i];
2006 if (delta > maxdelta)
2008 maxdelta = delta;
2011 // If any displacement is larger than the stepsize limit, reduce the step
2012 if (maxdelta > inputrec->em_stepsize)
2014 stepsize *= 0.1;
2017 while (maxdelta > inputrec->em_stepsize);
2019 // Take a trial step and move the coordinate array xc[] to position C
2020 for (i = 0; i < n; i++)
2022 xc[i] = lastx[i] + c*s[i];
2025 neval++;
2026 // Calculate energy for the trial step in position C
2027 ems.s.x = (rvec *)xc;
2028 ems.f = (rvec *)fc;
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 EpotC = ems.epot;
2036 // Calc line gradient in position C
2037 for (gpc = 0, i = 0; i < n; i++)
2039 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
2041 /* Sum the gradient along the line across CPUs */
2042 if (PAR(cr))
2044 gmx_sumd(1, &gpc, cr);
2047 // This is the max amount of increase in energy we tolerate.
2048 // By allowing VERY small changes (close to numerical precision) we
2049 // frequently find even better (lower) final energies.
2050 tmp = sqrt(GMX_REAL_EPS)*fabs(EpotA);
2052 // Accept the step if the energy is lower in the new position C (compared to A),
2053 // or if it is not significantly higher and the line derivative is still negative.
2054 if (EpotC < EpotA || (gpc < 0 && EpotC < (EpotA+tmp)))
2056 // Great, we found a better energy. We no longer try to alter the
2057 // stepsize, but simply accept this new better position. The we select a new
2058 // search direction instead, which will be much more efficient than continuing
2059 // to take smaller steps along a line. Set fnorm based on the new C position,
2060 // which will be used to update the stepsize to 1/fnorm further down.
2061 foundlower = TRUE;
2062 fnorm = ems.fnorm;
2064 else
2066 // If we got here, the energy is NOT lower in point C, i.e. it will be the same
2067 // or higher than in point A. In this case it is pointless to move to point C,
2068 // so we will have to do more iterations along the same line to find a smaller
2069 // value in the interval [A=0.0,C].
2070 // Here, A is still 0.0, but that will change when we do a search in the interval
2071 // [0.0,C] below. That search we will do by interpolation or bisection rather
2072 // than with the stepsize, so no need to modify it. For the next search direction
2073 // it will be reset to 1/fnorm anyway.
2074 foundlower = FALSE;
2077 if (!foundlower)
2079 // OK, if we didn't find a lower value we will have to locate one now - there must
2080 // be one in the interval [a,c].
2081 // The same thing is valid here, though: Don't spend dozens of iterations to find
2082 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2083 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2084 // I also have a safeguard for potentially really pathological functions so we never
2085 // take more than 20 steps before we give up.
2086 // If we already found a lower value we just skip this step and continue to the update.
2087 nminstep = 0;
2090 // Select a new trial point B in the interval [A,C].
2091 // If the derivatives at points a & c have different sign we interpolate to zero,
2092 // otherwise just do a bisection since there might be multiple minima/maxima
2093 // inside the interval.
2094 if (gpa < 0 && gpc > 0)
2096 b = a + gpa*(a-c)/(gpc-gpa);
2098 else
2100 b = 0.5*(a+c);
2103 /* safeguard if interpolation close to machine accuracy causes errors:
2104 * never go outside the interval
2106 if (b <= a || b >= c)
2108 b = 0.5*(a+c);
2111 // Take a trial step to point B
2112 for (i = 0; i < n; i++)
2114 xb[i] = lastx[i] + b*s[i];
2117 neval++;
2118 // Calculate energy for the trial step in point B
2119 ems.s.x = (rvec *)xb;
2120 ems.f = (rvec *)fb;
2121 evaluate_energy(fplog, cr,
2122 top_global, &ems, top,
2123 inputrec, nrnb, wcycle, gstat,
2124 vsite, constr, fcd, graph, mdatoms, fr,
2125 mu_tot, enerd, vir, pres, step, FALSE);
2126 EpotB = ems.epot;
2127 fnorm = ems.fnorm;
2129 // Calculate gradient in point B
2130 for (gpb = 0, i = 0; i < n; i++)
2132 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2135 /* Sum the gradient along the line across CPUs */
2136 if (PAR(cr))
2138 gmx_sumd(1, &gpb, cr);
2141 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2142 // at the new point B, and rename the endpoints of this new interval A and C.
2143 if (gpb > 0)
2145 /* Replace c endpoint with b */
2146 EpotC = EpotB;
2147 c = b;
2148 gpc = gpb;
2149 /* swap coord pointers b/c */
2150 xtmp = xb;
2151 ftmp = fb;
2152 xb = xc;
2153 fb = fc;
2154 xc = xtmp;
2155 fc = ftmp;
2157 else
2159 /* Replace a endpoint with b */
2160 EpotA = EpotB;
2161 a = b;
2162 gpa = gpb;
2163 /* swap coord pointers a/b */
2164 xtmp = xb;
2165 ftmp = fb;
2166 xb = xa;
2167 fb = fa;
2168 xa = xtmp;
2169 fa = ftmp;
2173 * Stop search as soon as we find a value smaller than the endpoints,
2174 * or if the tolerance is below machine precision.
2175 * Never run more than 20 steps, no matter what.
2177 nminstep++;
2179 while ((EpotB > EpotA || EpotB > EpotC) && (nminstep < 20));
2181 if (fabs(EpotB-Epot0) < GMX_REAL_EPS || nminstep >= 20)
2183 /* OK. We couldn't find a significantly lower energy.
2184 * If ncorr==0 this was steepest descent, and then we give up.
2185 * If not, reset memory to restart as steepest descent before quitting.
2187 if (ncorr == 0)
2189 /* Converged */
2190 converged = TRUE;
2191 break;
2193 else
2195 /* Reset memory */
2196 ncorr = 0;
2197 /* Search in gradient direction */
2198 for (i = 0; i < n; i++)
2200 dx[point][i] = ff[i];
2202 /* Reset stepsize */
2203 stepsize = 1.0/fnorm;
2204 continue;
2208 /* Select min energy state of A & C, put the best in xx/ff/Epot
2210 if (EpotC < EpotA)
2212 Epot = EpotC;
2213 /* Use state C */
2214 for (i = 0; i < n; i++)
2216 xx[i] = xc[i];
2217 ff[i] = fc[i];
2219 step_taken = c;
2221 else
2223 Epot = EpotA;
2224 /* Use state A */
2225 for (i = 0; i < n; i++)
2227 xx[i] = xa[i];
2228 ff[i] = fa[i];
2230 step_taken = a;
2234 else
2236 /* found lower */
2237 Epot = EpotC;
2238 /* Use state C */
2239 for (i = 0; i < n; i++)
2241 xx[i] = xc[i];
2242 ff[i] = fc[i];
2244 step_taken = c;
2247 /* Update the memory information, and calculate a new
2248 * approximation of the inverse hessian
2251 /* Have new data in Epot, xx, ff */
2252 if (ncorr < nmaxcorr)
2254 ncorr++;
2257 for (i = 0; i < n; i++)
2259 dg[point][i] = lastf[i]-ff[i];
2260 dx[point][i] *= step_taken;
2263 dgdg = 0;
2264 dgdx = 0;
2265 for (i = 0; i < n; i++)
2267 dgdg += dg[point][i]*dg[point][i];
2268 dgdx += dg[point][i]*dx[point][i];
2271 diag = dgdx/dgdg;
2273 rho[point] = 1.0/dgdx;
2274 point++;
2276 if (point >= nmaxcorr)
2278 point = 0;
2281 /* Update */
2282 for (i = 0; i < n; i++)
2284 p[i] = ff[i];
2287 cp = point;
2289 /* Recursive update. First go back over the memory points */
2290 for (k = 0; k < ncorr; k++)
2292 cp--;
2293 if (cp < 0)
2295 cp = ncorr-1;
2298 sq = 0;
2299 for (i = 0; i < n; i++)
2301 sq += dx[cp][i]*p[i];
2304 alpha[cp] = rho[cp]*sq;
2306 for (i = 0; i < n; i++)
2308 p[i] -= alpha[cp]*dg[cp][i];
2312 for (i = 0; i < n; i++)
2314 p[i] *= diag;
2317 /* And then go forward again */
2318 for (k = 0; k < ncorr; k++)
2320 yr = 0;
2321 for (i = 0; i < n; i++)
2323 yr += p[i]*dg[cp][i];
2326 beta = rho[cp]*yr;
2327 beta = alpha[cp]-beta;
2329 for (i = 0; i < n; i++)
2331 p[i] += beta*dx[cp][i];
2334 cp++;
2335 if (cp >= ncorr)
2337 cp = 0;
2341 for (i = 0; i < n; i++)
2343 if (!frozen[i])
2345 dx[point][i] = p[i];
2347 else
2349 dx[point][i] = 0;
2353 /* Test whether the convergence criterion is met */
2354 get_f_norm_max(cr, &(inputrec->opts), mdatoms, f, &fnorm, &fmax, &nfmax);
2356 /* Print it if necessary */
2357 if (MASTER(cr))
2359 if (bVerbose)
2361 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2362 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2363 step, Epot, fnorm/sqrtNumAtoms, fmax, nfmax+1);
2364 fflush(stderr);
2366 /* Store the new (lower) energies */
2367 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2368 mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
2369 NULL, NULL, vir, pres, NULL, mu_tot, constr);
2370 do_log = do_per_step(step, inputrec->nstlog);
2371 do_ene = do_per_step(step, inputrec->nstenergy);
2372 if (do_log)
2374 print_ebin_header(fplog, step, step);
2376 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2377 do_log ? fplog : NULL, step, step, eprNORMAL,
2378 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2381 /* Send x and E to IMD client, if bIMD is TRUE. */
2382 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
2384 IMD_send_positions(inputrec->imd);
2387 // Reset stepsize in we are doing more iterations
2388 stepsize = 1.0/fnorm;
2390 /* Stop when the maximum force lies below tolerance.
2391 * If we have reached machine precision, converged is already set to true.
2393 converged = converged || (fmax < inputrec->em_tol);
2395 } /* End of the loop */
2397 /* IMD cleanup, if bIMD is TRUE. */
2398 IMD_finalize(inputrec->bIMD, inputrec->imd);
2400 if (converged)
2402 step--; /* we never took that last step in this case */
2405 if (fmax > inputrec->em_tol)
2407 if (MASTER(cr))
2409 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2410 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2412 converged = FALSE;
2415 /* If we printed energy and/or logfile last step (which was the last step)
2416 * we don't have to do it again, but otherwise print the final values.
2418 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2420 print_ebin_header(fplog, step, step);
2422 if (!do_ene || !do_log) /* Write final energy file entries */
2424 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2425 !do_log ? fplog : NULL, step, step, eprNORMAL,
2426 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2429 /* Print some stuff... */
2430 if (MASTER(cr))
2432 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2435 /* IMPORTANT!
2436 * For accurate normal mode calculation it is imperative that we
2437 * store the last conformation into the full precision binary trajectory.
2439 * However, we should only do it if we did NOT already write this step
2440 * above (which we did if do_x or do_f was true).
2442 do_x = !do_per_step(step, inputrec->nstxout);
2443 do_f = !do_per_step(step, inputrec->nstfout);
2444 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2445 top_global, inputrec, step,
2446 &ems, state_global);
2448 if (MASTER(cr))
2450 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2451 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2452 number_steps, Epot, fmax, nfmax, fnorm/sqrtNumAtoms);
2453 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2454 number_steps, Epot, fmax, nfmax, fnorm/sqrtNumAtoms);
2456 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2459 finish_em(cr, outf, walltime_accounting, wcycle);
2461 /* To print the actual number of steps we needed somewhere */
2462 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2464 return 0;
2465 } /* That's all folks */
2467 /*! \brief Do steepest descents minimization
2468 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
2469 int nfile, const t_filenm fnm[],
2470 const gmx_output_env_t *oenv, gmx_bool bVerbose,
2471 int nstglobalcomm,
2472 gmx_vsite_t *vsite, gmx_constr_t constr,
2473 int stepout,
2474 t_inputrec *inputrec,
2475 gmx_mtop_t *top_global, t_fcdata *fcd,
2476 t_state *state_global,
2477 t_mdatoms *mdatoms,
2478 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2479 gmx_edsam_t ed,
2480 t_forcerec *fr,
2481 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
2482 real cpt_period, real max_hours,
2483 int imdport,
2484 unsigned long Flags,
2485 gmx_walltime_accounting_t walltime_accounting)
2487 double do_steep(FILE *fplog, t_commrec *cr, const gmx::MDLogger gmx_unused &mdlog,
2488 int nfile, const t_filenm fnm[],
2489 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
2490 int gmx_unused nstglobalcomm,
2491 gmx_vsite_t *vsite, gmx_constr_t constr,
2492 int gmx_unused stepout,
2493 t_inputrec *inputrec,
2494 gmx_mtop_t *top_global, t_fcdata *fcd,
2495 t_state *state_global,
2496 t_mdatoms *mdatoms,
2497 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2498 gmx_edsam_t gmx_unused ed,
2499 t_forcerec *fr,
2500 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2501 gmx_membed_t gmx_unused *membed,
2502 real gmx_unused cpt_period, real gmx_unused max_hours,
2503 int imdport,
2504 unsigned long gmx_unused Flags,
2505 gmx_walltime_accounting_t walltime_accounting)
2507 const char *SD = "Steepest Descents";
2508 em_state_t *s_min, *s_try;
2509 gmx_localtop_t *top;
2510 gmx_enerdata_t *enerd;
2511 rvec *f;
2512 gmx_global_stat_t gstat;
2513 t_graph *graph;
2514 real stepsize;
2515 real ustep, fnormn;
2516 gmx_mdoutf_t outf;
2517 t_mdebin *mdebin;
2518 gmx_bool bDone, bAbort, do_x, do_f;
2519 tensor vir, pres;
2520 rvec mu_tot;
2521 int nsteps;
2522 int count = 0;
2523 int steps_accepted = 0;
2525 s_min = init_em_state();
2526 s_try = init_em_state();
2528 /* Init em and store the local state in s_try */
2529 init_em(fplog, SD, cr, inputrec,
2530 state_global, top_global, s_try, &top, &f,
2531 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat,
2532 vsite, constr, NULL,
2533 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
2535 /* Print to log file */
2536 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2538 /* Set variables for stepsize (in nm). This is the largest
2539 * step that we are going to make in any direction.
2541 ustep = inputrec->em_stepsize;
2542 stepsize = 0;
2544 /* Max number of steps */
2545 nsteps = inputrec->nsteps;
2547 if (MASTER(cr))
2549 /* Print to the screen */
2550 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2552 if (fplog)
2554 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2557 /**** HERE STARTS THE LOOP ****
2558 * count is the counter for the number of steps
2559 * bDone will be TRUE when the minimization has converged
2560 * bAbort will be TRUE when nsteps steps have been performed or when
2561 * the stepsize becomes smaller than is reasonable for machine precision
2563 count = 0;
2564 bDone = FALSE;
2565 bAbort = FALSE;
2566 while (!bDone && !bAbort)
2568 bAbort = (nsteps >= 0) && (count == nsteps);
2570 /* set new coordinates, except for first step */
2571 bool validStep = true;
2572 if (count > 0)
2574 validStep =
2575 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2576 s_min, stepsize, s_min->f, s_try,
2577 constr, top, nrnb, wcycle, count);
2580 if (validStep)
2582 evaluate_energy(fplog, cr,
2583 top_global, s_try, top,
2584 inputrec, nrnb, wcycle, gstat,
2585 vsite, constr, fcd, graph, mdatoms, fr,
2586 mu_tot, enerd, vir, pres, count, count == 0);
2588 else
2590 // Signal constraint error during stepping with energy=inf
2591 s_try->epot = std::numeric_limits<real>::infinity();
2594 if (MASTER(cr))
2596 print_ebin_header(fplog, count, count);
2599 if (count == 0)
2601 s_min->epot = s_try->epot;
2604 /* Print it if necessary */
2605 if (MASTER(cr))
2607 if (bVerbose)
2609 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2610 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2611 ( (count == 0) || (s_try->epot < s_min->epot) ) ? '\n' : '\r');
2612 fflush(stderr);
2615 if ( (count == 0) || (s_try->epot < s_min->epot) )
2617 /* Store the new (lower) energies */
2618 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2619 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2620 s_try->s.box, NULL, NULL, vir, pres, NULL, mu_tot, constr);
2622 /* Prepare IMD energy record, if bIMD is TRUE. */
2623 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2625 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2626 do_per_step(steps_accepted, inputrec->nstdisreout),
2627 do_per_step(steps_accepted, inputrec->nstorireout),
2628 fplog, count, count, eprNORMAL,
2629 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2630 fflush(fplog);
2634 /* Now if the new energy is smaller than the previous...
2635 * or if this is the first step!
2636 * or if we did random steps!
2639 if ( (count == 0) || (s_try->epot < s_min->epot) )
2641 steps_accepted++;
2643 /* Test whether the convergence criterion is met... */
2644 bDone = (s_try->fmax < inputrec->em_tol);
2646 /* Copy the arrays for force, positions and energy */
2647 /* The 'Min' array always holds the coords and forces of the minimal
2648 sampled energy */
2649 swap_em_state(s_min, s_try);
2650 if (count > 0)
2652 ustep *= 1.2;
2655 /* Write to trn, if necessary */
2656 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2657 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2658 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
2659 top_global, inputrec, count,
2660 s_min, state_global);
2662 else
2664 /* If energy is not smaller make the step smaller... */
2665 ustep *= 0.5;
2667 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2669 /* Reload the old state */
2670 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2671 s_min, top, mdatoms, fr, vsite, constr,
2672 nrnb, wcycle);
2676 /* Determine new step */
2677 stepsize = ustep/s_min->fmax;
2679 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2680 #if GMX_DOUBLE
2681 if (count == nsteps || ustep < 1e-12)
2682 #else
2683 if (count == nsteps || ustep < 1e-6)
2684 #endif
2686 if (MASTER(cr))
2688 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2689 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2691 bAbort = TRUE;
2694 /* Send IMD energies and positions, if bIMD is TRUE. */
2695 if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
2697 IMD_send_positions(inputrec->imd);
2700 count++;
2701 } /* End of the loop */
2703 /* IMD cleanup, if bIMD is TRUE. */
2704 IMD_finalize(inputrec->bIMD, inputrec->imd);
2706 /* Print some data... */
2707 if (MASTER(cr))
2709 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2711 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2712 top_global, inputrec, count,
2713 s_min, state_global);
2715 if (MASTER(cr))
2717 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2718 fnormn = s_min->fnorm/sqrtNumAtoms;
2720 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2721 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2722 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2723 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2726 finish_em(cr, outf, walltime_accounting, wcycle);
2728 /* To print the actual number of steps we needed somewhere */
2729 inputrec->nsteps = count;
2731 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2733 return 0;
2734 } /* That's all folks */
2736 /*! \brief Do normal modes analysis
2737 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
2738 int nfile, const t_filenm fnm[],
2739 const gmx_output_env_t *oenv, gmx_bool bVerbose,
2740 int nstglobalcomm,
2741 gmx_vsite_t *vsite, gmx_constr_t constr,
2742 int stepout,
2743 t_inputrec *inputrec,
2744 gmx_mtop_t *top_global, t_fcdata *fcd,
2745 t_state *state_global,
2746 t_mdatoms *mdatoms,
2747 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2748 gmx_edsam_t ed,
2749 t_forcerec *fr,
2750 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
2751 real cpt_period, real max_hours,
2752 int imdport,
2753 unsigned long Flags,
2754 gmx_walltime_accounting_t walltime_accounting)
2756 double do_nm(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
2757 int nfile, const t_filenm fnm[],
2758 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
2759 int gmx_unused nstglobalcomm,
2760 gmx_vsite_t *vsite, gmx_constr_t constr,
2761 int gmx_unused stepout,
2762 t_inputrec *inputrec,
2763 gmx_mtop_t *top_global, t_fcdata *fcd,
2764 t_state *state_global,
2765 t_mdatoms *mdatoms,
2766 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2767 gmx_edsam_t gmx_unused ed,
2768 t_forcerec *fr,
2769 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2770 gmx_membed_t gmx_unused *membed,
2771 real gmx_unused cpt_period, real gmx_unused max_hours,
2772 int imdport,
2773 unsigned long gmx_unused Flags,
2774 gmx_walltime_accounting_t walltime_accounting)
2776 const char *NM = "Normal Mode Analysis";
2777 gmx_mdoutf_t outf;
2778 int nnodes, node;
2779 gmx_localtop_t *top;
2780 gmx_enerdata_t *enerd;
2781 rvec *f;
2782 gmx_global_stat_t gstat;
2783 t_graph *graph;
2784 tensor vir, pres;
2785 rvec mu_tot;
2786 rvec *fneg, *dfdx;
2787 gmx_bool bSparse; /* use sparse matrix storage format */
2788 size_t sz;
2789 gmx_sparsematrix_t * sparse_matrix = NULL;
2790 real * full_matrix = NULL;
2791 em_state_t * state_work;
2793 /* added with respect to mdrun */
2794 int row, col;
2795 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2796 real x_min;
2797 bool bIsMaster = MASTER(cr);
2799 if (constr != NULL)
2801 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2804 state_work = init_em_state();
2806 gmx_shellfc_t *shellfc;
2808 /* Init em and store the local state in state_minimum */
2809 init_em(fplog, NM, cr, inputrec,
2810 state_global, top_global, state_work, &top,
2812 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat,
2813 vsite, constr, &shellfc,
2814 nfile, fnm, &outf, NULL, imdport, Flags, wcycle);
2816 std::vector<size_t> atom_index = get_atom_index(top_global);
2817 snew(fneg, atom_index.size());
2818 snew(dfdx, atom_index.size());
2820 #if !GMX_DOUBLE
2821 if (bIsMaster)
2823 fprintf(stderr,
2824 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2825 " which MIGHT not be accurate enough for normal mode analysis.\n"
2826 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
2827 " are fairly modest even if you recompile in double precision.\n\n");
2829 #endif
2831 /* Check if we can/should use sparse storage format.
2833 * Sparse format is only useful when the Hessian itself is sparse, which it
2834 * will be when we use a cutoff.
2835 * For small systems (n<1000) it is easier to always use full matrix format, though.
2837 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2839 GMX_LOG(mdlog.warning).appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2840 bSparse = FALSE;
2842 else if (atom_index.size() < 1000)
2844 GMX_LOG(mdlog.warning).appendTextFormatted("Small system size (N=%d), using full Hessian format.",
2845 atom_index.size());
2846 bSparse = FALSE;
2848 else
2850 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2851 bSparse = TRUE;
2854 /* Number of dimensions, based on real atoms, that is not vsites or shell */
2855 sz = DIM*atom_index.size();
2857 fprintf(stderr, "Allocating Hessian memory...\n\n");
2859 if (bSparse)
2861 sparse_matrix = gmx_sparsematrix_init(sz);
2862 sparse_matrix->compressed_symmetric = TRUE;
2864 else
2866 snew(full_matrix, sz*sz);
2869 init_nrnb(nrnb);
2871 where();
2873 /* Write start time and temperature */
2874 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2876 /* fudge nr of steps to nr of atoms */
2877 inputrec->nsteps = atom_index.size()*2;
2879 if (bIsMaster)
2881 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2882 *(top_global->name), (int)inputrec->nsteps);
2885 nnodes = cr->nnodes;
2887 /* Make evaluate_energy do a single node force calculation */
2888 cr->nnodes = 1;
2889 evaluate_energy(fplog, cr,
2890 top_global, state_work, top,
2891 inputrec, nrnb, wcycle, gstat,
2892 vsite, constr, fcd, graph, mdatoms, fr,
2893 mu_tot, enerd, vir, pres, -1, TRUE);
2894 cr->nnodes = nnodes;
2896 /* if forces are not small, warn user */
2897 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, state_work);
2899 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work->fmax);
2900 if (state_work->fmax > 1.0e-3)
2902 GMX_LOG(mdlog.warning).appendText(
2903 "The force is probably not small enough to "
2904 "ensure that you are at a minimum.\n"
2905 "Be aware that negative eigenvalues may occur\n"
2906 "when the resulting matrix is diagonalized.");
2909 /***********************************************************
2911 * Loop over all pairs in matrix
2913 * do_force called twice. Once with positive and
2914 * once with negative displacement
2916 ************************************************************/
2918 /* Steps are divided one by one over the nodes */
2919 bool bNS = true;
2920 for (unsigned int aid = cr->nodeid; aid < atom_index.size(); aid += nnodes)
2922 size_t atom = atom_index[aid];
2923 for (size_t d = 0; d < DIM; d++)
2925 gmx_bool bBornRadii = FALSE;
2926 gmx_int64_t step = 0;
2927 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2928 double t = 0;
2930 x_min = state_work->s.x[atom][d];
2932 for (unsigned int dx = 0; (dx < 2); dx++)
2934 if (dx == 0)
2936 state_work->s.x[atom][d] = x_min - der_range;
2938 else
2940 state_work->s.x[atom][d] = x_min + der_range;
2943 /* Make evaluate_energy do a single node force calculation */
2944 cr->nnodes = 1;
2945 if (shellfc)
2947 /* Now is the time to relax the shells */
2948 (void) relax_shell_flexcon(fplog, cr, bVerbose, step,
2949 inputrec, bNS, force_flags,
2950 top,
2951 constr, enerd, fcd,
2952 &state_work->s, state_work->f, vir, mdatoms,
2953 nrnb, wcycle, graph, &top_global->groups,
2954 shellfc, fr, bBornRadii, t, mu_tot,
2955 vsite, NULL);
2956 bNS = false;
2957 step++;
2959 else
2961 evaluate_energy(fplog, cr,
2962 top_global, state_work, top,
2963 inputrec, nrnb, wcycle, gstat,
2964 vsite, constr, fcd, graph, mdatoms, fr,
2965 mu_tot, enerd, vir, pres, atom*2+dx, FALSE);
2968 cr->nnodes = nnodes;
2970 if (dx == 0)
2972 for (size_t i = 0; i < atom_index.size(); i++)
2974 copy_rvec(state_work->f[atom_index[i]], fneg[i]);
2979 /* x is restored to original */
2980 state_work->s.x[atom][d] = x_min;
2982 for (size_t j = 0; j < atom_index.size(); j++)
2984 for (size_t k = 0; (k < DIM); k++)
2986 dfdx[j][k] =
2987 -(state_work->f[atom_index[j]][k] - fneg[j][k])/(2*der_range);
2991 if (!bIsMaster)
2993 #if GMX_MPI
2994 #define mpi_type GMX_MPI_REAL
2995 MPI_Send(dfdx[0], atom_index.size()*DIM, mpi_type, MASTER(cr),
2996 cr->nodeid, cr->mpi_comm_mygroup);
2997 #endif
2999 else
3001 for (node = 0; (node < nnodes && atom+node < atom_index.size()); node++)
3003 if (node > 0)
3005 #if GMX_MPI
3006 MPI_Status stat;
3007 MPI_Recv(dfdx[0], atom_index.size()*DIM, mpi_type, node, node,
3008 cr->mpi_comm_mygroup, &stat);
3009 #undef mpi_type
3010 #endif
3013 row = (atom + node)*DIM + d;
3015 for (size_t j = 0; j < atom_index.size(); j++)
3017 for (size_t k = 0; k < DIM; k++)
3019 col = j*DIM + k;
3021 if (bSparse)
3023 if (col >= row && dfdx[j][k] != 0.0)
3025 gmx_sparsematrix_increment_value(sparse_matrix,
3026 row, col, dfdx[j][k]);
3029 else
3031 full_matrix[row*sz+col] = dfdx[j][k];
3038 if (bVerbose && fplog)
3040 fflush(fplog);
3043 /* write progress */
3044 if (bIsMaster && bVerbose)
3046 fprintf(stderr, "\rFinished step %d out of %d",
3047 static_cast<int>(std::min(atom+nnodes, atom_index.size())),
3048 static_cast<int>(atom_index.size()));
3049 fflush(stderr);
3053 if (bIsMaster)
3055 fprintf(stderr, "\n\nWriting Hessian...\n");
3056 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
3059 finish_em(cr, outf, walltime_accounting, wcycle);
3061 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size()*2);
3063 return 0;
3066 } // namespace gmx