Prepared t_mdatoms for using vector
[gromacs.git] / src / gromacs / mdlib / minimize.cpp
blob16c8b2c8659aedb0725174b46c9d7e3571bb0d89
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,2017, 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/mdtypes/state.h"
91 #include "gromacs/pbcutil/mshift.h"
92 #include "gromacs/pbcutil/pbc.h"
93 #include "gromacs/timing/wallcycle.h"
94 #include "gromacs/timing/walltime_accounting.h"
95 #include "gromacs/topology/mtop_util.h"
96 #include "gromacs/topology/topology.h"
97 #include "gromacs/utility/cstringutil.h"
98 #include "gromacs/utility/exceptions.h"
99 #include "gromacs/utility/fatalerror.h"
100 #include "gromacs/utility/logger.h"
101 #include "gromacs/utility/smalloc.h"
103 //! Utility structure for manipulating states during EM
104 typedef struct {
105 //! Copy of the global state
106 t_state s;
107 //! Force array
108 PaddedRVecVector f;
109 //! Potential energy
110 real epot;
111 //! Norm of the force
112 real fnorm;
113 //! Maximum force
114 real fmax;
115 //! Direction
116 int a_fmax;
117 } em_state_t;
119 //! Print the EM starting conditions
120 static void print_em_start(FILE *fplog,
121 t_commrec *cr,
122 gmx_walltime_accounting_t walltime_accounting,
123 gmx_wallcycle_t wcycle,
124 const char *name)
126 walltime_accounting_start(walltime_accounting);
127 wallcycle_start(wcycle, ewcRUN);
128 print_start(fplog, cr, walltime_accounting, name);
131 //! Stop counting time for EM
132 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
133 gmx_wallcycle_t wcycle)
135 wallcycle_stop(wcycle, ewcRUN);
137 walltime_accounting_end(walltime_accounting);
140 //! Printing a log file and console header
141 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
143 fprintf(out, "\n");
144 fprintf(out, "%s:\n", minimizer);
145 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
146 fprintf(out, " Number of steps = %12d\n", nsteps);
149 //! Print warning message
150 static void warn_step(FILE *fp, real ftol, gmx_bool bLastStep, gmx_bool bConstrain)
152 char buffer[2048];
153 if (bLastStep)
155 sprintf(buffer,
156 "\nEnergy minimization reached the maximum number "
157 "of steps before the forces reached the requested "
158 "precision Fmax < %g.\n", ftol);
160 else
162 sprintf(buffer,
163 "\nEnergy minimization has stopped, but the forces have "
164 "not converged to the requested precision Fmax < %g (which "
165 "may not be possible for your system). It stopped "
166 "because the algorithm tried to make a new step whose size "
167 "was too small, or there was no change in the energy since "
168 "last step. Either way, we regard the minimization as "
169 "converged to within the available machine precision, "
170 "given your starting configuration and EM parameters.\n%s%s",
171 ftol,
172 sizeof(real) < sizeof(double) ?
173 "\nDouble precision normally gives you higher accuracy, but "
174 "this is often not needed for preparing to run molecular "
175 "dynamics.\n" :
177 bConstrain ?
178 "You might need to increase your constraint accuracy, or turn\n"
179 "off constraints altogether (set constraints = none in mdp file)\n" :
180 "");
182 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
185 //! Print message about convergence of the EM
186 static void print_converged(FILE *fp, const char *alg, real ftol,
187 gmx_int64_t count, gmx_bool bDone, gmx_int64_t nsteps,
188 const em_state_t *ems, double sqrtNumAtoms)
190 char buf[STEPSTRSIZE];
192 if (bDone)
194 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
195 alg, ftol, gmx_step_str(count, buf));
197 else if (count < nsteps)
199 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
200 "but did not reach the requested Fmax < %g.\n",
201 alg, gmx_step_str(count, buf), ftol);
203 else
205 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
206 alg, ftol, gmx_step_str(count, buf));
209 #if GMX_DOUBLE
210 fprintf(fp, "Potential Energy = %21.14e\n", ems->epot);
211 fprintf(fp, "Maximum force = %21.14e on atom %d\n", ems->fmax, ems->a_fmax + 1);
212 fprintf(fp, "Norm of force = %21.14e\n", ems->fnorm/sqrtNumAtoms);
213 #else
214 fprintf(fp, "Potential Energy = %14.7e\n", ems->epot);
215 fprintf(fp, "Maximum force = %14.7e on atom %d\n", ems->fmax, ems->a_fmax + 1);
216 fprintf(fp, "Norm of force = %14.7e\n", ems->fnorm/sqrtNumAtoms);
217 #endif
220 //! Compute the norm and max of the force array in parallel
221 static void get_f_norm_max(t_commrec *cr,
222 t_grpopts *opts, t_mdatoms *mdatoms, const rvec *f,
223 real *fnorm, real *fmax, int *a_fmax)
225 double fnorm2, *sum;
226 real fmax2, fam;
227 int la_max, a_max, start, end, i, m, gf;
229 /* This routine finds the largest force and returns it.
230 * On parallel machines the global max is taken.
232 fnorm2 = 0;
233 fmax2 = 0;
234 la_max = -1;
235 start = 0;
236 end = mdatoms->homenr;
237 if (mdatoms->cFREEZE)
239 for (i = start; i < end; i++)
241 gf = mdatoms->cFREEZE[i];
242 fam = 0;
243 for (m = 0; m < DIM; m++)
245 if (!opts->nFreeze[gf][m])
247 fam += gmx::square(f[i][m]);
250 fnorm2 += fam;
251 if (fam > fmax2)
253 fmax2 = fam;
254 la_max = i;
258 else
260 for (i = start; i < end; i++)
262 fam = norm2(f[i]);
263 fnorm2 += fam;
264 if (fam > fmax2)
266 fmax2 = fam;
267 la_max = i;
272 if (la_max >= 0 && DOMAINDECOMP(cr))
274 a_max = cr->dd->gatindex[la_max];
276 else
278 a_max = la_max;
280 if (PAR(cr))
282 snew(sum, 2*cr->nnodes+1);
283 sum[2*cr->nodeid] = fmax2;
284 sum[2*cr->nodeid+1] = a_max;
285 sum[2*cr->nnodes] = fnorm2;
286 gmx_sumd(2*cr->nnodes+1, sum, cr);
287 fnorm2 = sum[2*cr->nnodes];
288 /* Determine the global maximum */
289 for (i = 0; i < cr->nnodes; i++)
291 if (sum[2*i] > fmax2)
293 fmax2 = sum[2*i];
294 a_max = (int)(sum[2*i+1] + 0.5);
297 sfree(sum);
300 if (fnorm)
302 *fnorm = sqrt(fnorm2);
304 if (fmax)
306 *fmax = sqrt(fmax2);
308 if (a_fmax)
310 *a_fmax = a_max;
314 //! Compute the norm of the force
315 static void get_state_f_norm_max(t_commrec *cr,
316 t_grpopts *opts, t_mdatoms *mdatoms,
317 em_state_t *ems)
319 get_f_norm_max(cr, opts, mdatoms, as_rvec_array(ems->f.data()),
320 &ems->fnorm, &ems->fmax, &ems->a_fmax);
323 //! Initialize the energy minimization
324 static void init_em(FILE *fplog, const char *title,
325 t_commrec *cr, gmx::IMDOutputProvider *outputProvider,
326 t_inputrec *ir,
327 const MdrunOptions &mdrunOptions,
328 t_state *state_global, gmx_mtop_t *top_global,
329 em_state_t *ems, gmx_localtop_t **top,
330 t_nrnb *nrnb, rvec mu_tot,
331 t_forcerec *fr, gmx_enerdata_t **enerd,
332 t_graph **graph, gmx::MDAtoms *mdAtoms, gmx_global_stat_t *gstat,
333 gmx_vsite_t *vsite, gmx_constr_t constr, gmx_shellfc_t **shellfc,
334 int nfile, const t_filenm fnm[],
335 gmx_mdoutf_t *outf, t_mdebin **mdebin,
336 gmx_wallcycle_t wcycle)
338 real dvdl_constr;
340 if (fplog)
342 fprintf(fplog, "Initiating %s\n", title);
345 if (MASTER(cr))
347 state_global->ngtc = 0;
349 /* Initialize lambda variables */
350 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, nullptr);
353 init_nrnb(nrnb);
355 /* Interactive molecular dynamics */
356 init_IMD(ir, cr, top_global, fplog, 1,
357 MASTER(cr) ? as_rvec_array(state_global->x.data()) : nullptr,
358 nfile, fnm, nullptr, mdrunOptions);
360 if (ir->eI == eiNM)
362 GMX_ASSERT(shellfc != NULL, "With NM we always support shells");
364 *shellfc = init_shell_flexcon(stdout,
365 top_global,
366 n_flexible_constraints(constr),
367 ir->nstcalcenergy,
368 DOMAINDECOMP(cr));
370 else
372 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI), "This else currently only handles energy minimizers, consider if your algorithm needs shell/flexible-constraint support");
374 /* With energy minimization, shells and flexible constraints are
375 * automatically minimized when treated like normal DOFS.
377 if (shellfc != nullptr)
379 *shellfc = nullptr;
383 auto mdatoms = mdAtoms->mdatoms();
384 if (DOMAINDECOMP(cr))
386 *top = dd_init_local_top(top_global);
388 dd_init_local_state(cr->dd, state_global, &ems->s);
390 /* Distribute the charge groups over the nodes from the master node */
391 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
392 state_global, top_global, ir,
393 &ems->s, &ems->f, mdAtoms, *top,
394 fr, vsite, constr,
395 nrnb, nullptr, FALSE);
396 dd_store_state(cr->dd, &ems->s);
398 *graph = nullptr;
400 else
402 state_change_natoms(state_global, state_global->natoms);
403 /* Just copy the state */
404 ems->s = *state_global;
405 state_change_natoms(&ems->s, ems->s.natoms);
406 /* We need to allocate one element extra, since we might use
407 * (unaligned) 4-wide SIMD loads to access rvec entries.
409 ems->f.resize(gmx::paddedRVecVectorSize(ems->s.natoms));
411 snew(*top, 1);
412 mdAlgorithmsSetupAtomData(cr, ir, top_global, *top, fr,
413 graph, mdAtoms,
414 vsite, shellfc ? *shellfc : nullptr);
416 if (vsite)
418 set_vsite_top(vsite, *top, mdatoms);
422 update_mdatoms(mdAtoms->mdatoms(), ems->s.lambda[efptMASS]);
424 if (constr)
426 if (ir->eConstrAlg == econtSHAKE &&
427 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
429 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
430 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
433 if (!DOMAINDECOMP(cr))
435 set_constraints(constr, *top, ir, mdatoms, cr);
438 if (!ir->bContinuation)
440 /* Constrain the starting coordinates */
441 dvdl_constr = 0;
442 constrain(PAR(cr) ? nullptr : fplog, TRUE, TRUE, constr, &(*top)->idef,
443 ir, cr, -1, 0, 1.0, mdatoms,
444 as_rvec_array(ems->s.x.data()),
445 as_rvec_array(ems->s.x.data()),
446 nullptr,
447 fr->bMolPBC, ems->s.box,
448 ems->s.lambda[efptFEP], &dvdl_constr,
449 nullptr, nullptr, nrnb, econqCoord);
453 if (PAR(cr))
455 *gstat = global_stat_init(ir);
457 else
459 *gstat = nullptr;
462 *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, ir, top_global, nullptr, wcycle);
464 snew(*enerd, 1);
465 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
466 *enerd);
468 if (mdebin != nullptr)
470 /* Init bin for energy stuff */
471 *mdebin = init_mdebin(mdoutf_get_fp_ene(*outf), top_global, ir, nullptr);
474 clear_rvec(mu_tot);
475 calc_shifts(ems->s.box, fr->shift_vec);
478 //! Finalize the minimization
479 static void finish_em(t_commrec *cr, gmx_mdoutf_t outf,
480 gmx_walltime_accounting_t walltime_accounting,
481 gmx_wallcycle_t wcycle)
483 if (!thisRankHasDuty(cr, DUTY_PME))
485 /* Tell the PME only node to finish */
486 gmx_pme_send_finish(cr);
489 done_mdoutf(outf);
491 em_time_end(walltime_accounting, wcycle);
494 //! Swap two different EM states during minimization
495 static void swap_em_state(em_state_t **ems1, em_state_t **ems2)
497 em_state_t *tmp;
499 tmp = *ems1;
500 *ems1 = *ems2;
501 *ems2 = tmp;
504 //! Save the EM trajectory
505 static void write_em_traj(FILE *fplog, t_commrec *cr,
506 gmx_mdoutf_t outf,
507 gmx_bool bX, gmx_bool bF, const char *confout,
508 gmx_mtop_t *top_global,
509 t_inputrec *ir, gmx_int64_t step,
510 em_state_t *state,
511 t_state *state_global,
512 ObservablesHistory *observablesHistory)
514 int mdof_flags = 0;
516 if (bX)
518 mdof_flags |= MDOF_X;
520 if (bF)
522 mdof_flags |= MDOF_F;
525 /* If we want IMD output, set appropriate MDOF flag */
526 if (ir->bIMD)
528 mdof_flags |= MDOF_IMD;
531 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
532 top_global, step, (double)step,
533 &state->s, state_global, observablesHistory,
534 &state->f);
536 if (confout != nullptr && MASTER(cr))
538 GMX_RELEASE_ASSERT(bX, "The code below assumes that (with domain decomposition), x is collected to state_global in the call above.");
539 /* With domain decomposition the call above collected the state->s.x
540 * into state_global->x. Without DD we copy the local state pointer.
542 if (!DOMAINDECOMP(cr))
544 state_global = &state->s;
547 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
549 /* Make molecules whole only for confout writing */
550 do_pbc_mtop(fplog, ir->ePBC, state->s.box, top_global,
551 as_rvec_array(state_global->x.data()));
554 write_sto_conf_mtop(confout,
555 *top_global->name, top_global,
556 as_rvec_array(state_global->x.data()), nullptr, ir->ePBC, state->s.box);
560 //! \brief Do one minimization step
562 // \returns true when the step succeeded, false when a constraint error occurred
563 static bool do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
564 gmx_bool bMolPBC,
565 em_state_t *ems1, real a, const PaddedRVecVector *force,
566 em_state_t *ems2,
567 gmx_constr_t constr, gmx_localtop_t *top,
568 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
569 gmx_int64_t count)
572 t_state *s1, *s2;
573 int start, end;
574 real dvdl_constr;
575 int nthreads gmx_unused;
577 bool validStep = true;
579 s1 = &ems1->s;
580 s2 = &ems2->s;
582 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
584 gmx_incons("state mismatch in do_em_step");
587 s2->flags = s1->flags;
589 if (s2->natoms != s1->natoms)
591 state_change_natoms(s2, s1->natoms);
592 /* We need to allocate one element extra, since we might use
593 * (unaligned) 4-wide SIMD loads to access rvec entries.
595 ems2->f.resize(gmx::paddedRVecVectorSize(s2->natoms));
597 if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
599 s2->cg_gl.resize(s1->cg_gl.size());
602 copy_mat(s1->box, s2->box);
603 /* Copy free energy state */
604 s2->lambda = s1->lambda;
605 copy_mat(s1->box, s2->box);
607 start = 0;
608 end = md->homenr;
610 // cppcheck-suppress unreadVariable
611 nthreads = gmx_omp_nthreads_get(emntUpdate);
612 #pragma omp parallel num_threads(nthreads)
614 const rvec *x1 = as_rvec_array(s1->x.data());
615 rvec *x2 = as_rvec_array(s2->x.data());
616 const rvec *f = as_rvec_array(force->data());
618 int gf = 0;
619 #pragma omp for schedule(static) nowait
620 for (int i = start; i < end; i++)
624 if (md->cFREEZE)
626 gf = md->cFREEZE[i];
628 for (int m = 0; m < DIM; m++)
630 if (ir->opts.nFreeze[gf][m])
632 x2[i][m] = x1[i][m];
634 else
636 x2[i][m] = x1[i][m] + a*f[i][m];
640 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
643 if (s2->flags & (1<<estCGP))
645 /* Copy the CG p vector */
646 const rvec *p1 = as_rvec_array(s1->cg_p.data());
647 rvec *p2 = as_rvec_array(s2->cg_p.data());
648 #pragma omp for schedule(static) nowait
649 for (int i = start; i < end; i++)
651 // Trivial OpenMP block that does not throw
652 copy_rvec(p1[i], p2[i]);
656 if (DOMAINDECOMP(cr))
658 s2->ddp_count = s1->ddp_count;
660 /* OpenMP does not supported unsigned loop variables */
661 #pragma omp for schedule(static) nowait
662 for (int i = 0; i < static_cast<int>(s2->cg_gl.size()); i++)
664 s2->cg_gl[i] = s1->cg_gl[i];
666 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
670 if (constr)
672 wallcycle_start(wcycle, ewcCONSTR);
673 dvdl_constr = 0;
674 validStep =
675 constrain(nullptr, TRUE, TRUE, constr, &top->idef,
676 ir, cr, count, 0, 1.0, md,
677 as_rvec_array(s1->x.data()), as_rvec_array(s2->x.data()),
678 nullptr, bMolPBC, s2->box,
679 s2->lambda[efptBONDED], &dvdl_constr,
680 nullptr, nullptr, nrnb, econqCoord);
681 wallcycle_stop(wcycle, ewcCONSTR);
683 // We should move this check to the different minimizers
684 if (!validStep && ir->eI != eiSteep)
686 gmx_fatal(FARGS, "The coordinates could not be constrained. Minimizer '%s' can not handle constraint failures, use minimizer '%s' before using '%s'.",
687 EI(ir->eI), EI(eiSteep), EI(ir->eI));
691 return validStep;
694 //! Prepare EM for using domain decomposition parallellization
695 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
696 gmx_mtop_t *top_global, t_inputrec *ir,
697 em_state_t *ems, gmx_localtop_t *top,
698 gmx::MDAtoms *mdAtoms, t_forcerec *fr,
699 gmx_vsite_t *vsite, gmx_constr_t constr,
700 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
702 /* Repartition the domain decomposition */
703 dd_partition_system(fplog, step, cr, FALSE, 1,
704 nullptr, top_global, ir,
705 &ems->s, &ems->f,
706 mdAtoms, top, fr, vsite, constr,
707 nrnb, wcycle, FALSE);
708 dd_store_state(cr->dd, &ems->s);
711 //! De one energy evaluation
712 static void evaluate_energy(FILE *fplog, t_commrec *cr,
713 gmx_mtop_t *top_global,
714 em_state_t *ems, gmx_localtop_t *top,
715 t_inputrec *inputrec,
716 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
717 gmx_global_stat_t gstat,
718 gmx_vsite_t *vsite, gmx_constr_t constr,
719 t_fcdata *fcd,
720 t_graph *graph, gmx::MDAtoms *mdAtoms,
721 t_forcerec *fr, rvec mu_tot,
722 gmx_enerdata_t *enerd, tensor vir, tensor pres,
723 gmx_int64_t count, gmx_bool bFirst)
725 real t;
726 gmx_bool bNS;
727 tensor force_vir, shake_vir, ekin;
728 real dvdl_constr, prescorr, enercorr, dvdlcorr;
729 real terminate = 0;
731 /* Set the time to the initial time, the time does not change during EM */
732 t = inputrec->init_t;
734 if (bFirst ||
735 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
737 /* This is the first state or an old state used before the last ns */
738 bNS = TRUE;
740 else
742 bNS = FALSE;
743 if (inputrec->nstlist > 0)
745 bNS = TRUE;
749 if (vsite)
751 construct_vsites(vsite, as_rvec_array(ems->s.x.data()), 1, nullptr,
752 top->idef.iparams, top->idef.il,
753 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
756 if (DOMAINDECOMP(cr) && bNS)
758 /* Repartition the domain decomposition */
759 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
760 ems, top, mdAtoms, fr, vsite, constr,
761 nrnb, wcycle);
764 /* Calc force & energy on new trial position */
765 /* do_force always puts the charge groups in the box and shifts again
766 * We do not unshift, so molecules are always whole in congrad.c
768 do_force(fplog, cr, inputrec,
769 count, nrnb, wcycle, top, &top_global->groups,
770 ems->s.box, &ems->s.x, &ems->s.hist,
771 &ems->f, force_vir, mdAtoms->mdatoms(), enerd, fcd,
772 ems->s.lambda, graph, fr, vsite, mu_tot, t, nullptr, TRUE,
773 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
774 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
775 (bNS ? GMX_FORCE_NS : 0),
776 DOMAINDECOMP(cr) ?
777 DdOpenBalanceRegionBeforeForceComputation::yes :
778 DdOpenBalanceRegionBeforeForceComputation::no,
779 DOMAINDECOMP(cr) ?
780 DdCloseBalanceRegionAfterForceComputation::yes :
781 DdCloseBalanceRegionAfterForceComputation::no);
783 /* Clear the unused shake virial and pressure */
784 clear_mat(shake_vir);
785 clear_mat(pres);
787 /* Communicate stuff when parallel */
788 if (PAR(cr) && inputrec->eI != eiNM)
790 wallcycle_start(wcycle, ewcMoveE);
792 global_stat(gstat, cr, enerd, force_vir, shake_vir, mu_tot,
793 inputrec, nullptr, nullptr, nullptr, 1, &terminate,
794 nullptr, FALSE,
795 CGLO_ENERGY |
796 CGLO_PRESSURE |
797 CGLO_CONSTRAINT);
799 wallcycle_stop(wcycle, ewcMoveE);
802 /* Calculate long range corrections to pressure and energy */
803 calc_dispcorr(inputrec, fr, ems->s.box, ems->s.lambda[efptVDW],
804 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
805 enerd->term[F_DISPCORR] = enercorr;
806 enerd->term[F_EPOT] += enercorr;
807 enerd->term[F_PRES] += prescorr;
808 enerd->term[F_DVDL] += dvdlcorr;
810 ems->epot = enerd->term[F_EPOT];
812 if (constr)
814 /* Project out the constraint components of the force */
815 wallcycle_start(wcycle, ewcCONSTR);
816 dvdl_constr = 0;
817 rvec *f_rvec = as_rvec_array(ems->f.data());
818 constrain(nullptr, FALSE, FALSE, constr, &top->idef,
819 inputrec, cr, count, 0, 1.0, mdAtoms->mdatoms(),
820 as_rvec_array(ems->s.x.data()), f_rvec, f_rvec,
821 fr->bMolPBC, ems->s.box,
822 ems->s.lambda[efptBONDED], &dvdl_constr,
823 nullptr, &shake_vir, nrnb, econqForceDispl);
824 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
825 m_add(force_vir, shake_vir, vir);
826 wallcycle_stop(wcycle, ewcCONSTR);
828 else
830 copy_mat(force_vir, vir);
833 clear_mat(ekin);
834 enerd->term[F_PRES] =
835 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
837 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
839 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
841 get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
845 //! Parallel utility summing energies and forces
846 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
847 gmx_mtop_t *top_global,
848 em_state_t *s_min, em_state_t *s_b)
850 t_block *cgs_gl;
851 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
852 double partsum;
853 unsigned char *grpnrFREEZE;
855 if (debug)
857 fprintf(debug, "Doing reorder_partsum\n");
860 const rvec *fm = as_rvec_array(s_min->f.data());
861 const rvec *fb = as_rvec_array(s_b->f.data());
863 cgs_gl = dd_charge_groups_global(cr->dd);
864 index = cgs_gl->index;
866 /* Collect fm in a global vector fmg.
867 * This conflicts with the spirit of domain decomposition,
868 * but to fully optimize this a much more complicated algorithm is required.
870 rvec *fmg;
871 snew(fmg, top_global->natoms);
873 ncg = s_min->s.cg_gl.size();
874 cg_gl = s_min->s.cg_gl.data();
875 i = 0;
876 for (c = 0; c < ncg; c++)
878 cg = cg_gl[c];
879 a0 = index[cg];
880 a1 = index[cg+1];
881 for (a = a0; a < a1; a++)
883 copy_rvec(fm[i], fmg[a]);
884 i++;
887 gmx_sum(top_global->natoms*3, fmg[0], cr);
889 /* Now we will determine the part of the sum for the cgs in state s_b */
890 ncg = s_b->s.cg_gl.size();
891 cg_gl = s_b->s.cg_gl.data();
892 partsum = 0;
893 i = 0;
894 gf = 0;
895 grpnrFREEZE = top_global->groups.grpnr[egcFREEZE];
896 for (c = 0; c < ncg; c++)
898 cg = cg_gl[c];
899 a0 = index[cg];
900 a1 = index[cg+1];
901 for (a = a0; a < a1; a++)
903 if (mdatoms->cFREEZE && grpnrFREEZE)
905 gf = grpnrFREEZE[i];
907 for (m = 0; m < DIM; m++)
909 if (!opts->nFreeze[gf][m])
911 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
914 i++;
918 sfree(fmg);
920 return partsum;
923 //! Print some stuff, like beta, whatever that means.
924 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
925 gmx_mtop_t *top_global,
926 em_state_t *s_min, em_state_t *s_b)
928 double sum;
930 /* This is just the classical Polak-Ribiere calculation of beta;
931 * it looks a bit complicated since we take freeze groups into account,
932 * and might have to sum it in parallel runs.
935 if (!DOMAINDECOMP(cr) ||
936 (s_min->s.ddp_count == cr->dd->ddp_count &&
937 s_b->s.ddp_count == cr->dd->ddp_count))
939 const rvec *fm = as_rvec_array(s_min->f.data());
940 const rvec *fb = as_rvec_array(s_b->f.data());
941 sum = 0;
942 int gf = 0;
943 /* This part of code can be incorrect with DD,
944 * since the atom ordering in s_b and s_min might differ.
946 for (int i = 0; i < mdatoms->homenr; i++)
948 if (mdatoms->cFREEZE)
950 gf = mdatoms->cFREEZE[i];
952 for (int m = 0; m < DIM; m++)
954 if (!opts->nFreeze[gf][m])
956 sum += (fb[i][m] - fm[i][m])*fb[i][m];
961 else
963 /* We need to reorder cgs while summing */
964 sum = reorder_partsum(cr, opts, mdatoms, top_global, s_min, s_b);
966 if (PAR(cr))
968 gmx_sumd(1, &sum, cr);
971 return sum/gmx::square(s_min->fnorm);
974 namespace gmx
977 /*! \brief Do conjugate gradients minimization
978 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
979 int nfile, const t_filenm fnm[],
980 const gmx_output_env_t *oenv,
981 const MdrunOptions &mdrunOptions,
982 gmx_vsite_t *vsite, gmx_constr_t constr,
983 gmx::IMDOutputProvider *outputProvider,
984 t_inputrec *inputrec,
985 gmx_mtop_t *top_global, t_fcdata *fcd,
986 t_state *state_global,
987 gmx::MDAtoms *mdAtoms,
988 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
989 gmx_edsam_t ed,
990 t_forcerec *fr,
991 const ReplicaExchangeParameters &replExParams,
992 gmx_membed_t gmx_unused *membed,
993 gmx_walltime_accounting_t walltime_accounting)
995 double do_cg(FILE *fplog, t_commrec *cr, const gmx::MDLogger gmx_unused &mdlog,
996 int nfile, const t_filenm fnm[],
997 const gmx_output_env_t gmx_unused *oenv,
998 const MdrunOptions &mdrunOptions,
999 gmx_vsite_t *vsite, gmx_constr_t constr,
1000 gmx::IMDOutputProvider *outputProvider,
1001 t_inputrec *inputrec,
1002 gmx_mtop_t *top_global, t_fcdata *fcd,
1003 t_state *state_global,
1004 ObservablesHistory *observablesHistory,
1005 gmx::MDAtoms *mdAtoms,
1006 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1007 t_forcerec *fr,
1008 const ReplicaExchangeParameters gmx_unused &replExParams,
1009 gmx_membed_t gmx_unused *membed,
1010 gmx_walltime_accounting_t walltime_accounting)
1012 const char *CG = "Polak-Ribiere Conjugate Gradients";
1014 gmx_localtop_t *top;
1015 gmx_enerdata_t *enerd;
1016 gmx_global_stat_t gstat;
1017 t_graph *graph;
1018 double tmp, minstep;
1019 real stepsize;
1020 real a, b, c, beta = 0.0;
1021 real epot_repl = 0;
1022 real pnorm;
1023 t_mdebin *mdebin;
1024 gmx_bool converged, foundlower;
1025 rvec mu_tot;
1026 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1027 tensor vir, pres;
1028 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1029 gmx_mdoutf_t outf;
1030 int m, step, nminstep;
1031 auto mdatoms = mdAtoms->mdatoms();
1033 step = 0;
1035 // Ensure the extra per-atom state array gets allocated
1036 state_global->flags |= (1<<estCGP);
1038 /* Create 4 states on the stack and extract pointers that we will swap */
1039 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1040 em_state_t *s_min = &s0;
1041 em_state_t *s_a = &s1;
1042 em_state_t *s_b = &s2;
1043 em_state_t *s_c = &s3;
1045 /* Init em and store the local state in s_min */
1046 init_em(fplog, CG, cr, outputProvider, inputrec, mdrunOptions,
1047 state_global, top_global, s_min, &top,
1048 nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
1049 vsite, constr, nullptr,
1050 nfile, fnm, &outf, &mdebin, wcycle);
1052 /* Print to log file */
1053 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1055 /* Max number of steps */
1056 number_steps = inputrec->nsteps;
1058 if (MASTER(cr))
1060 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1062 if (fplog)
1064 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1067 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1068 /* do_force always puts the charge groups in the box and shifts again
1069 * We do not unshift, so molecules are always whole in congrad.c
1071 evaluate_energy(fplog, cr,
1072 top_global, s_min, top,
1073 inputrec, nrnb, wcycle, gstat,
1074 vsite, constr, fcd, graph, mdAtoms, fr,
1075 mu_tot, enerd, vir, pres, -1, TRUE);
1076 where();
1078 if (MASTER(cr))
1080 /* Copy stuff to the energy bin for easy printing etc. */
1081 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1082 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1083 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1085 print_ebin_header(fplog, step, step);
1086 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1087 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1089 where();
1091 /* Estimate/guess the initial stepsize */
1092 stepsize = inputrec->em_stepsize/s_min->fnorm;
1094 if (MASTER(cr))
1096 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1097 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1098 s_min->fmax, s_min->a_fmax+1);
1099 fprintf(stderr, " F-Norm = %12.5e\n",
1100 s_min->fnorm/sqrtNumAtoms);
1101 fprintf(stderr, "\n");
1102 /* and copy to the log file too... */
1103 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1104 s_min->fmax, s_min->a_fmax+1);
1105 fprintf(fplog, " F-Norm = %12.5e\n",
1106 s_min->fnorm/sqrtNumAtoms);
1107 fprintf(fplog, "\n");
1109 /* Start the loop over CG steps.
1110 * Each successful step is counted, and we continue until
1111 * we either converge or reach the max number of steps.
1113 converged = FALSE;
1114 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1117 /* start taking steps in a new direction
1118 * First time we enter the routine, beta=0, and the direction is
1119 * simply the negative gradient.
1122 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1123 rvec *pm = as_rvec_array(s_min->s.cg_p.data());
1124 const rvec *sfm = as_rvec_array(s_min->f.data());
1125 double gpa = 0;
1126 int gf = 0;
1127 for (int i = 0; i < mdatoms->homenr; i++)
1129 if (mdatoms->cFREEZE)
1131 gf = mdatoms->cFREEZE[i];
1133 for (m = 0; m < DIM; m++)
1135 if (!inputrec->opts.nFreeze[gf][m])
1137 pm[i][m] = sfm[i][m] + beta*pm[i][m];
1138 gpa -= pm[i][m]*sfm[i][m];
1139 /* f is negative gradient, thus the sign */
1141 else
1143 pm[i][m] = 0;
1148 /* Sum the gradient along the line across CPUs */
1149 if (PAR(cr))
1151 gmx_sumd(1, &gpa, cr);
1154 /* Calculate the norm of the search vector */
1155 get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1157 /* Just in case stepsize reaches zero due to numerical precision... */
1158 if (stepsize <= 0)
1160 stepsize = inputrec->em_stepsize/pnorm;
1164 * Double check the value of the derivative in the search direction.
1165 * If it is positive it must be due to the old information in the
1166 * CG formula, so just remove that and start over with beta=0.
1167 * This corresponds to a steepest descent step.
1169 if (gpa > 0)
1171 beta = 0;
1172 step--; /* Don't count this step since we are restarting */
1173 continue; /* Go back to the beginning of the big for-loop */
1176 /* Calculate minimum allowed stepsize, before the average (norm)
1177 * relative change in coordinate is smaller than precision
1179 minstep = 0;
1180 for (int i = 0; i < mdatoms->homenr; i++)
1182 for (m = 0; m < DIM; m++)
1184 tmp = fabs(s_min->s.x[i][m]);
1185 if (tmp < 1.0)
1187 tmp = 1.0;
1189 tmp = pm[i][m]/tmp;
1190 minstep += tmp*tmp;
1193 /* Add up from all CPUs */
1194 if (PAR(cr))
1196 gmx_sumd(1, &minstep, cr);
1199 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1201 if (stepsize < minstep)
1203 converged = TRUE;
1204 break;
1207 /* Write coordinates if necessary */
1208 do_x = do_per_step(step, inputrec->nstxout);
1209 do_f = do_per_step(step, inputrec->nstfout);
1211 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
1212 top_global, inputrec, step,
1213 s_min, state_global, observablesHistory);
1215 /* Take a step downhill.
1216 * In theory, we should minimize the function along this direction.
1217 * That is quite possible, but it turns out to take 5-10 function evaluations
1218 * for each line. However, we dont really need to find the exact minimum -
1219 * it is much better to start a new CG step in a modified direction as soon
1220 * as we are close to it. This will save a lot of energy evaluations.
1222 * In practice, we just try to take a single step.
1223 * If it worked (i.e. lowered the energy), we increase the stepsize but
1224 * the continue straight to the next CG step without trying to find any minimum.
1225 * If it didn't work (higher energy), there must be a minimum somewhere between
1226 * the old position and the new one.
1228 * Due to the finite numerical accuracy, it turns out that it is a good idea
1229 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1230 * This leads to lower final energies in the tests I've done. / Erik
1232 s_a->epot = s_min->epot;
1233 a = 0.0;
1234 c = a + stepsize; /* reference position along line is zero */
1236 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1238 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1239 s_min, top, mdAtoms, fr, vsite, constr,
1240 nrnb, wcycle);
1243 /* Take a trial step (new coords in s_c) */
1244 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, &s_min->s.cg_p, s_c,
1245 constr, top, nrnb, wcycle, -1);
1247 neval++;
1248 /* Calculate energy for the trial step */
1249 evaluate_energy(fplog, cr,
1250 top_global, s_c, top,
1251 inputrec, nrnb, wcycle, gstat,
1252 vsite, constr, fcd, graph, mdAtoms, fr,
1253 mu_tot, enerd, vir, pres, -1, FALSE);
1255 /* Calc derivative along line */
1256 const rvec *pc = as_rvec_array(s_c->s.cg_p.data());
1257 const rvec *sfc = as_rvec_array(s_c->f.data());
1258 double gpc = 0;
1259 for (int i = 0; i < mdatoms->homenr; i++)
1261 for (m = 0; m < DIM; m++)
1263 gpc -= pc[i][m]*sfc[i][m]; /* f is negative gradient, thus the sign */
1266 /* Sum the gradient along the line across CPUs */
1267 if (PAR(cr))
1269 gmx_sumd(1, &gpc, cr);
1272 /* This is the max amount of increase in energy we tolerate */
1273 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1275 /* Accept the step if the energy is lower, or if it is not significantly higher
1276 * and the line derivative is still negative.
1278 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1280 foundlower = TRUE;
1281 /* Great, we found a better energy. Increase step for next iteration
1282 * if we are still going down, decrease it otherwise
1284 if (gpc < 0)
1286 stepsize *= 1.618034; /* The golden section */
1288 else
1290 stepsize *= 0.618034; /* 1/golden section */
1293 else
1295 /* New energy is the same or higher. We will have to do some work
1296 * to find a smaller value in the interval. Take smaller step next time!
1298 foundlower = FALSE;
1299 stepsize *= 0.618034;
1305 /* OK, if we didn't find a lower value we will have to locate one now - there must
1306 * be one in the interval [a=0,c].
1307 * The same thing is valid here, though: Don't spend dozens of iterations to find
1308 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1309 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1311 * I also have a safeguard for potentially really pathological functions so we never
1312 * take more than 20 steps before we give up ...
1314 * If we already found a lower value we just skip this step and continue to the update.
1316 double gpb;
1317 if (!foundlower)
1319 nminstep = 0;
1323 /* Select a new trial point.
1324 * If the derivatives at points a & c have different sign we interpolate to zero,
1325 * otherwise just do a bisection.
1327 if (gpa < 0 && gpc > 0)
1329 b = a + gpa*(a-c)/(gpc-gpa);
1331 else
1333 b = 0.5*(a+c);
1336 /* safeguard if interpolation close to machine accuracy causes errors:
1337 * never go outside the interval
1339 if (b <= a || b >= c)
1341 b = 0.5*(a+c);
1344 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1346 /* Reload the old state */
1347 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1348 s_min, top, mdAtoms, fr, vsite, constr,
1349 nrnb, wcycle);
1352 /* Take a trial step to this new point - new coords in s_b */
1353 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, &s_min->s.cg_p, s_b,
1354 constr, top, nrnb, wcycle, -1);
1356 neval++;
1357 /* Calculate energy for the trial step */
1358 evaluate_energy(fplog, cr,
1359 top_global, s_b, top,
1360 inputrec, nrnb, wcycle, gstat,
1361 vsite, constr, fcd, graph, mdAtoms, fr,
1362 mu_tot, enerd, vir, pres, -1, FALSE);
1364 /* p does not change within a step, but since the domain decomposition
1365 * might change, we have to use cg_p of s_b here.
1367 const rvec *pb = as_rvec_array(s_b->s.cg_p.data());
1368 const rvec *sfb = as_rvec_array(s_b->f.data());
1369 gpb = 0;
1370 for (int i = 0; i < mdatoms->homenr; i++)
1372 for (m = 0; m < DIM; m++)
1374 gpb -= pb[i][m]*sfb[i][m]; /* f is negative gradient, thus the sign */
1377 /* Sum the gradient along the line across CPUs */
1378 if (PAR(cr))
1380 gmx_sumd(1, &gpb, cr);
1383 if (debug)
1385 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1386 s_a->epot, s_b->epot, s_c->epot, gpb);
1389 epot_repl = s_b->epot;
1391 /* Keep one of the intervals based on the value of the derivative at the new point */
1392 if (gpb > 0)
1394 /* Replace c endpoint with b */
1395 swap_em_state(&s_b, &s_c);
1396 c = b;
1397 gpc = gpb;
1399 else
1401 /* Replace a endpoint with b */
1402 swap_em_state(&s_b, &s_a);
1403 a = b;
1404 gpa = gpb;
1408 * Stop search as soon as we find a value smaller than the endpoints.
1409 * Never run more than 20 steps, no matter what.
1411 nminstep++;
1413 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1414 (nminstep < 20));
1416 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1417 nminstep >= 20)
1419 /* OK. We couldn't find a significantly lower energy.
1420 * If beta==0 this was steepest descent, and then we give up.
1421 * If not, set beta=0 and restart with steepest descent before quitting.
1423 if (beta == 0.0)
1425 /* Converged */
1426 converged = TRUE;
1427 break;
1429 else
1431 /* Reset memory before giving up */
1432 beta = 0.0;
1433 continue;
1437 /* Select min energy state of A & C, put the best in B.
1439 if (s_c->epot < s_a->epot)
1441 if (debug)
1443 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1444 s_c->epot, s_a->epot);
1446 swap_em_state(&s_b, &s_c);
1447 gpb = gpc;
1449 else
1451 if (debug)
1453 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1454 s_a->epot, s_c->epot);
1456 swap_em_state(&s_b, &s_a);
1457 gpb = gpa;
1461 else
1463 if (debug)
1465 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1466 s_c->epot);
1468 swap_em_state(&s_b, &s_c);
1469 gpb = gpc;
1472 /* new search direction */
1473 /* beta = 0 means forget all memory and restart with steepest descents. */
1474 if (nstcg && ((step % nstcg) == 0))
1476 beta = 0.0;
1478 else
1480 /* s_min->fnorm cannot be zero, because then we would have converged
1481 * and broken out.
1484 /* Polak-Ribiere update.
1485 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1487 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1489 /* Limit beta to prevent oscillations */
1490 if (fabs(beta) > 5.0)
1492 beta = 0.0;
1496 /* update positions */
1497 swap_em_state(&s_min, &s_b);
1498 gpa = gpb;
1500 /* Print it if necessary */
1501 if (MASTER(cr))
1503 if (mdrunOptions.verbose)
1505 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1506 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1507 step, s_min->epot, s_min->fnorm/sqrtNumAtoms,
1508 s_min->fmax, s_min->a_fmax+1);
1509 fflush(stderr);
1511 /* Store the new (lower) energies */
1512 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1513 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1514 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1516 do_log = do_per_step(step, inputrec->nstlog);
1517 do_ene = do_per_step(step, inputrec->nstenergy);
1519 /* Prepare IMD energy record, if bIMD is TRUE. */
1520 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1522 if (do_log)
1524 print_ebin_header(fplog, step, step);
1526 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1527 do_log ? fplog : nullptr, step, step, eprNORMAL,
1528 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1531 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1532 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, as_rvec_array(state_global->x.data()), inputrec, 0, wcycle) && MASTER(cr))
1534 IMD_send_positions(inputrec->imd);
1537 /* Stop when the maximum force lies below tolerance.
1538 * If we have reached machine precision, converged is already set to true.
1540 converged = converged || (s_min->fmax < inputrec->em_tol);
1542 } /* End of the loop */
1544 /* IMD cleanup, if bIMD is TRUE. */
1545 IMD_finalize(inputrec->bIMD, inputrec->imd);
1547 if (converged)
1549 step--; /* we never took that last step in this case */
1552 if (s_min->fmax > inputrec->em_tol)
1554 if (MASTER(cr))
1556 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1557 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1559 converged = FALSE;
1562 if (MASTER(cr))
1564 /* If we printed energy and/or logfile last step (which was the last step)
1565 * we don't have to do it again, but otherwise print the final values.
1567 if (!do_log)
1569 /* Write final value to log since we didn't do anything the last step */
1570 print_ebin_header(fplog, step, step);
1572 if (!do_ene || !do_log)
1574 /* Write final energy file entries */
1575 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1576 !do_log ? fplog : nullptr, step, step, eprNORMAL,
1577 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1581 /* Print some stuff... */
1582 if (MASTER(cr))
1584 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1587 /* IMPORTANT!
1588 * For accurate normal mode calculation it is imperative that we
1589 * store the last conformation into the full precision binary trajectory.
1591 * However, we should only do it if we did NOT already write this step
1592 * above (which we did if do_x or do_f was true).
1594 do_x = !do_per_step(step, inputrec->nstxout);
1595 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1597 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1598 top_global, inputrec, step,
1599 s_min, state_global, observablesHistory);
1602 if (MASTER(cr))
1604 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1605 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1606 s_min, sqrtNumAtoms);
1607 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1608 s_min, sqrtNumAtoms);
1610 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1613 finish_em(cr, outf, walltime_accounting, wcycle);
1615 /* To print the actual number of steps we needed somewhere */
1616 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1618 return 0;
1619 } /* That's all folks */
1622 /*! \brief Do L-BFGS conjugate gradients minimization
1623 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
1624 int nfile, const t_filenm fnm[],
1625 const gmx_output_env_t *oenv,
1626 const MdrunOptions &mdrunOptions,
1627 gmx_vsite_t *vsite, gmx_constr_t constr,
1628 gmx::IMDOutputProvider *outputProvider,
1629 t_inputrec *inputrec,
1630 gmx_mtop_t *top_global, t_fcdata *fcd,
1631 t_state *state_global,
1632 gmx::MDAtoms *mdAtoms,
1633 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1634 gmx_edsam_t ed,
1635 t_forcerec *fr,
1636 const ReplicaExchangeParameters &replExParams,
1637 gmx_membed_t gmx_unused *membed,
1638 gmx_walltime_accounting_t walltime_accounting)
1640 double do_lbfgs(FILE *fplog, t_commrec *cr, const gmx::MDLogger gmx_unused &mdlog,
1641 int nfile, const t_filenm fnm[],
1642 const gmx_output_env_t gmx_unused *oenv,
1643 const MdrunOptions &mdrunOptions,
1644 gmx_vsite_t *vsite, gmx_constr_t constr,
1645 gmx::IMDOutputProvider *outputProvider,
1646 t_inputrec *inputrec,
1647 gmx_mtop_t *top_global, t_fcdata *fcd,
1648 t_state *state_global,
1649 ObservablesHistory *observablesHistory,
1650 gmx::MDAtoms *mdAtoms,
1651 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1652 t_forcerec *fr,
1653 const ReplicaExchangeParameters gmx_unused &replExParams,
1654 gmx_membed_t gmx_unused *membed,
1655 gmx_walltime_accounting_t walltime_accounting)
1657 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1658 em_state_t ems;
1659 gmx_localtop_t *top;
1660 gmx_enerdata_t *enerd;
1661 gmx_global_stat_t gstat;
1662 t_graph *graph;
1663 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1664 double stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1665 real *rho, *alpha, *p, *s, **dx, **dg;
1666 real a, b, c, maxdelta, delta;
1667 real diag, Epot0;
1668 real dgdx, dgdg, sq, yr, beta;
1669 t_mdebin *mdebin;
1670 gmx_bool converged;
1671 rvec mu_tot;
1672 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1673 tensor vir, pres;
1674 int start, end, number_steps;
1675 gmx_mdoutf_t outf;
1676 int i, k, m, n, gf, step;
1677 int mdof_flags;
1678 auto mdatoms = mdAtoms->mdatoms();
1680 if (PAR(cr))
1682 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1685 if (nullptr != constr)
1687 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).");
1690 n = 3*state_global->natoms;
1691 nmaxcorr = inputrec->nbfgscorr;
1693 snew(frozen, n);
1695 snew(p, n);
1696 snew(rho, nmaxcorr);
1697 snew(alpha, nmaxcorr);
1699 snew(dx, nmaxcorr);
1700 for (i = 0; i < nmaxcorr; i++)
1702 snew(dx[i], n);
1705 snew(dg, nmaxcorr);
1706 for (i = 0; i < nmaxcorr; i++)
1708 snew(dg[i], n);
1711 step = 0;
1712 neval = 0;
1714 /* Init em */
1715 init_em(fplog, LBFGS, cr, outputProvider, inputrec, mdrunOptions,
1716 state_global, top_global, &ems, &top,
1717 nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
1718 vsite, constr, nullptr,
1719 nfile, fnm, &outf, &mdebin, wcycle);
1721 start = 0;
1722 end = mdatoms->homenr;
1724 /* We need 4 working states */
1725 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1726 em_state_t *sa = &s0;
1727 em_state_t *sb = &s1;
1728 em_state_t *sc = &s2;
1729 em_state_t *last = &s3;
1730 /* Initialize by copying the state from ems (we could skip x and f here) */
1731 *sa = ems;
1732 *sb = ems;
1733 *sc = ems;
1735 /* Print to log file */
1736 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1738 do_log = do_ene = do_x = do_f = TRUE;
1740 /* Max number of steps */
1741 number_steps = inputrec->nsteps;
1743 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1744 gf = 0;
1745 for (i = start; i < end; i++)
1747 if (mdatoms->cFREEZE)
1749 gf = mdatoms->cFREEZE[i];
1751 for (m = 0; m < DIM; m++)
1753 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1756 if (MASTER(cr))
1758 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1760 if (fplog)
1762 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1765 if (vsite)
1767 construct_vsites(vsite, as_rvec_array(state_global->x.data()), 1, nullptr,
1768 top->idef.iparams, top->idef.il,
1769 fr->ePBC, fr->bMolPBC, cr, state_global->box);
1772 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1773 /* do_force always puts the charge groups in the box and shifts again
1774 * We do not unshift, so molecules are always whole
1776 neval++;
1777 evaluate_energy(fplog, cr,
1778 top_global, &ems, top,
1779 inputrec, nrnb, wcycle, gstat,
1780 vsite, constr, fcd, graph, mdAtoms, fr,
1781 mu_tot, enerd, vir, pres, -1, TRUE);
1782 where();
1784 if (MASTER(cr))
1786 /* Copy stuff to the energy bin for easy printing etc. */
1787 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1788 mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
1789 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1791 print_ebin_header(fplog, step, step);
1792 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1793 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1795 where();
1797 /* Set the initial step.
1798 * since it will be multiplied by the non-normalized search direction
1799 * vector (force vector the first time), we scale it by the
1800 * norm of the force.
1803 if (MASTER(cr))
1805 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1806 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1807 fprintf(stderr, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1808 fprintf(stderr, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1809 fprintf(stderr, "\n");
1810 /* and copy to the log file too... */
1811 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1812 fprintf(fplog, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1813 fprintf(fplog, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1814 fprintf(fplog, "\n");
1817 // Point is an index to the memory of search directions, where 0 is the first one.
1818 point = 0;
1820 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1821 real *fInit = static_cast<real *>(as_rvec_array(ems.f.data())[0]);
1822 for (i = 0; i < n; i++)
1824 if (!frozen[i])
1826 dx[point][i] = fInit[i]; /* Initial search direction */
1828 else
1830 dx[point][i] = 0;
1834 // Stepsize will be modified during the search, and actually it is not critical
1835 // (the main efficiency in the algorithm comes from changing directions), but
1836 // we still need an initial value, so estimate it as the inverse of the norm
1837 // so we take small steps where the potential fluctuates a lot.
1838 stepsize = 1.0/ems.fnorm;
1840 /* Start the loop over BFGS steps.
1841 * Each successful step is counted, and we continue until
1842 * we either converge or reach the max number of steps.
1845 ncorr = 0;
1847 /* Set the gradient from the force */
1848 converged = FALSE;
1849 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1852 /* Write coordinates if necessary */
1853 do_x = do_per_step(step, inputrec->nstxout);
1854 do_f = do_per_step(step, inputrec->nstfout);
1856 mdof_flags = 0;
1857 if (do_x)
1859 mdof_flags |= MDOF_X;
1862 if (do_f)
1864 mdof_flags |= MDOF_F;
1867 if (inputrec->bIMD)
1869 mdof_flags |= MDOF_IMD;
1872 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1873 top_global, step, (real)step, &ems.s, state_global, observablesHistory, &ems.f);
1875 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1877 /* make s a pointer to current search direction - point=0 first time we get here */
1878 s = dx[point];
1880 real *xx = static_cast<real *>(as_rvec_array(ems.s.x.data())[0]);
1881 real *ff = static_cast<real *>(as_rvec_array(ems.f.data())[0]);
1883 // calculate line gradient in position A
1884 for (gpa = 0, i = 0; i < n; i++)
1886 gpa -= s[i]*ff[i];
1889 /* Calculate minimum allowed stepsize along the line, before the average (norm)
1890 * relative change in coordinate is smaller than precision
1892 for (minstep = 0, i = 0; i < n; i++)
1894 tmp = fabs(xx[i]);
1895 if (tmp < 1.0)
1897 tmp = 1.0;
1899 tmp = s[i]/tmp;
1900 minstep += tmp*tmp;
1902 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1904 if (stepsize < minstep)
1906 converged = TRUE;
1907 break;
1910 // Before taking any steps along the line, store the old position
1911 *last = ems;
1912 real *lastx = static_cast<real *>(as_rvec_array(last->s.x.data())[0]);
1913 real *lastf = static_cast<real *>(as_rvec_array(last->f.data())[0]);
1914 Epot0 = ems.epot;
1916 *sa = ems;
1918 /* Take a step downhill.
1919 * In theory, we should find the actual minimum of the function in this
1920 * direction, somewhere along the line.
1921 * That is quite possible, but it turns out to take 5-10 function evaluations
1922 * for each line. However, we dont really need to find the exact minimum -
1923 * it is much better to start a new BFGS step in a modified direction as soon
1924 * as we are close to it. This will save a lot of energy evaluations.
1926 * In practice, we just try to take a single step.
1927 * If it worked (i.e. lowered the energy), we increase the stepsize but
1928 * continue straight to the next BFGS step without trying to find any minimum,
1929 * i.e. we change the search direction too. If the line was smooth, it is
1930 * likely we are in a smooth region, and then it makes sense to take longer
1931 * steps in the modified search direction too.
1933 * If it didn't work (higher energy), there must be a minimum somewhere between
1934 * the old position and the new one. Then we need to start by finding a lower
1935 * value before we change search direction. Since the energy was apparently
1936 * quite rough, we need to decrease the step size.
1938 * Due to the finite numerical accuracy, it turns out that it is a good idea
1939 * to accept a SMALL increase in energy, if the derivative is still downhill.
1940 * This leads to lower final energies in the tests I've done. / Erik
1943 // State "A" is the first position along the line.
1944 // reference position along line is initially zero
1945 a = 0.0;
1947 // Check stepsize first. We do not allow displacements
1948 // larger than emstep.
1952 // Pick a new position C by adding stepsize to A.
1953 c = a + stepsize;
1955 // Calculate what the largest change in any individual coordinate
1956 // would be (translation along line * gradient along line)
1957 maxdelta = 0;
1958 for (i = 0; i < n; i++)
1960 delta = c*s[i];
1961 if (delta > maxdelta)
1963 maxdelta = delta;
1966 // If any displacement is larger than the stepsize limit, reduce the step
1967 if (maxdelta > inputrec->em_stepsize)
1969 stepsize *= 0.1;
1972 while (maxdelta > inputrec->em_stepsize);
1974 // Take a trial step and move the coordinate array xc[] to position C
1975 real *xc = static_cast<real *>(as_rvec_array(sc->s.x.data())[0]);
1976 for (i = 0; i < n; i++)
1978 xc[i] = lastx[i] + c*s[i];
1981 neval++;
1982 // Calculate energy for the trial step in position C
1983 evaluate_energy(fplog, cr,
1984 top_global, sc, top,
1985 inputrec, nrnb, wcycle, gstat,
1986 vsite, constr, fcd, graph, mdAtoms, fr,
1987 mu_tot, enerd, vir, pres, step, FALSE);
1989 // Calc line gradient in position C
1990 real *fc = static_cast<real *>(as_rvec_array(sc->f.data())[0]);
1991 for (gpc = 0, i = 0; i < n; i++)
1993 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1995 /* Sum the gradient along the line across CPUs */
1996 if (PAR(cr))
1998 gmx_sumd(1, &gpc, cr);
2001 // This is the max amount of increase in energy we tolerate.
2002 // By allowing VERY small changes (close to numerical precision) we
2003 // frequently find even better (lower) final energies.
2004 tmp = sqrt(GMX_REAL_EPS)*fabs(sa->epot);
2006 // Accept the step if the energy is lower in the new position C (compared to A),
2007 // or if it is not significantly higher and the line derivative is still negative.
2008 if (sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp)))
2010 // Great, we found a better energy. We no longer try to alter the
2011 // stepsize, but simply accept this new better position. The we select a new
2012 // search direction instead, which will be much more efficient than continuing
2013 // to take smaller steps along a line. Set fnorm based on the new C position,
2014 // which will be used to update the stepsize to 1/fnorm further down.
2015 foundlower = TRUE;
2017 else
2019 // If we got here, the energy is NOT lower in point C, i.e. it will be the same
2020 // or higher than in point A. In this case it is pointless to move to point C,
2021 // so we will have to do more iterations along the same line to find a smaller
2022 // value in the interval [A=0.0,C].
2023 // Here, A is still 0.0, but that will change when we do a search in the interval
2024 // [0.0,C] below. That search we will do by interpolation or bisection rather
2025 // than with the stepsize, so no need to modify it. For the next search direction
2026 // it will be reset to 1/fnorm anyway.
2027 foundlower = FALSE;
2030 if (!foundlower)
2032 // OK, if we didn't find a lower value we will have to locate one now - there must
2033 // be one in the interval [a,c].
2034 // The same thing is valid here, though: Don't spend dozens of iterations to find
2035 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2036 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2037 // I also have a safeguard for potentially really pathological functions so we never
2038 // take more than 20 steps before we give up.
2039 // If we already found a lower value we just skip this step and continue to the update.
2040 real fnorm = 0;
2041 nminstep = 0;
2044 // Select a new trial point B in the interval [A,C].
2045 // If the derivatives at points a & c have different sign we interpolate to zero,
2046 // otherwise just do a bisection since there might be multiple minima/maxima
2047 // inside the interval.
2048 if (gpa < 0 && gpc > 0)
2050 b = a + gpa*(a-c)/(gpc-gpa);
2052 else
2054 b = 0.5*(a+c);
2057 /* safeguard if interpolation close to machine accuracy causes errors:
2058 * never go outside the interval
2060 if (b <= a || b >= c)
2062 b = 0.5*(a+c);
2065 // Take a trial step to point B
2066 real *xb = static_cast<real *>(as_rvec_array(sb->s.x.data())[0]);
2067 for (i = 0; i < n; i++)
2069 xb[i] = lastx[i] + b*s[i];
2072 neval++;
2073 // Calculate energy for the trial step in point B
2074 evaluate_energy(fplog, cr,
2075 top_global, sb, top,
2076 inputrec, nrnb, wcycle, gstat,
2077 vsite, constr, fcd, graph, mdAtoms, fr,
2078 mu_tot, enerd, vir, pres, step, FALSE);
2079 fnorm = sb->fnorm;
2081 // Calculate gradient in point B
2082 real *fb = static_cast<real *>(as_rvec_array(sb->f.data())[0]);
2083 for (gpb = 0, i = 0; i < n; i++)
2085 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2088 /* Sum the gradient along the line across CPUs */
2089 if (PAR(cr))
2091 gmx_sumd(1, &gpb, cr);
2094 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2095 // at the new point B, and rename the endpoints of this new interval A and C.
2096 if (gpb > 0)
2098 /* Replace c endpoint with b */
2099 c = b;
2100 /* swap states b and c */
2101 swap_em_state(&sb, &sc);
2103 else
2105 /* Replace a endpoint with b */
2106 a = b;
2107 /* swap states a and b */
2108 swap_em_state(&sa, &sb);
2112 * Stop search as soon as we find a value smaller than the endpoints,
2113 * or if the tolerance is below machine precision.
2114 * Never run more than 20 steps, no matter what.
2116 nminstep++;
2118 while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2120 if (fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2122 /* OK. We couldn't find a significantly lower energy.
2123 * If ncorr==0 this was steepest descent, and then we give up.
2124 * If not, reset memory to restart as steepest descent before quitting.
2126 if (ncorr == 0)
2128 /* Converged */
2129 converged = TRUE;
2130 break;
2132 else
2134 /* Reset memory */
2135 ncorr = 0;
2136 /* Search in gradient direction */
2137 for (i = 0; i < n; i++)
2139 dx[point][i] = ff[i];
2141 /* Reset stepsize */
2142 stepsize = 1.0/fnorm;
2143 continue;
2147 /* Select min energy state of A & C, put the best in xx/ff/Epot
2149 if (sc->epot < sa->epot)
2151 /* Use state C */
2152 ems = *sc;
2153 step_taken = c;
2155 else
2157 /* Use state A */
2158 ems = *sa;
2159 step_taken = a;
2163 else
2165 /* found lower */
2166 /* Use state C */
2167 ems = *sc;
2168 step_taken = c;
2171 /* Update the memory information, and calculate a new
2172 * approximation of the inverse hessian
2175 /* Have new data in Epot, xx, ff */
2176 if (ncorr < nmaxcorr)
2178 ncorr++;
2181 for (i = 0; i < n; i++)
2183 dg[point][i] = lastf[i]-ff[i];
2184 dx[point][i] *= step_taken;
2187 dgdg = 0;
2188 dgdx = 0;
2189 for (i = 0; i < n; i++)
2191 dgdg += dg[point][i]*dg[point][i];
2192 dgdx += dg[point][i]*dx[point][i];
2195 diag = dgdx/dgdg;
2197 rho[point] = 1.0/dgdx;
2198 point++;
2200 if (point >= nmaxcorr)
2202 point = 0;
2205 /* Update */
2206 for (i = 0; i < n; i++)
2208 p[i] = ff[i];
2211 cp = point;
2213 /* Recursive update. First go back over the memory points */
2214 for (k = 0; k < ncorr; k++)
2216 cp--;
2217 if (cp < 0)
2219 cp = ncorr-1;
2222 sq = 0;
2223 for (i = 0; i < n; i++)
2225 sq += dx[cp][i]*p[i];
2228 alpha[cp] = rho[cp]*sq;
2230 for (i = 0; i < n; i++)
2232 p[i] -= alpha[cp]*dg[cp][i];
2236 for (i = 0; i < n; i++)
2238 p[i] *= diag;
2241 /* And then go forward again */
2242 for (k = 0; k < ncorr; k++)
2244 yr = 0;
2245 for (i = 0; i < n; i++)
2247 yr += p[i]*dg[cp][i];
2250 beta = rho[cp]*yr;
2251 beta = alpha[cp]-beta;
2253 for (i = 0; i < n; i++)
2255 p[i] += beta*dx[cp][i];
2258 cp++;
2259 if (cp >= ncorr)
2261 cp = 0;
2265 for (i = 0; i < n; i++)
2267 if (!frozen[i])
2269 dx[point][i] = p[i];
2271 else
2273 dx[point][i] = 0;
2277 /* Print it if necessary */
2278 if (MASTER(cr))
2280 if (mdrunOptions.verbose)
2282 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2283 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2284 step, ems.epot, ems.fnorm/sqrtNumAtoms, ems.fmax, ems.a_fmax + 1);
2285 fflush(stderr);
2287 /* Store the new (lower) energies */
2288 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2289 mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
2290 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2291 do_log = do_per_step(step, inputrec->nstlog);
2292 do_ene = do_per_step(step, inputrec->nstenergy);
2293 if (do_log)
2295 print_ebin_header(fplog, step, step);
2297 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2298 do_log ? fplog : nullptr, step, step, eprNORMAL,
2299 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2302 /* Send x and E to IMD client, if bIMD is TRUE. */
2303 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, as_rvec_array(state_global->x.data()), inputrec, 0, wcycle) && MASTER(cr))
2305 IMD_send_positions(inputrec->imd);
2308 // Reset stepsize in we are doing more iterations
2309 stepsize = 1.0/ems.fnorm;
2311 /* Stop when the maximum force lies below tolerance.
2312 * If we have reached machine precision, converged is already set to true.
2314 converged = converged || (ems.fmax < inputrec->em_tol);
2316 } /* End of the loop */
2318 /* IMD cleanup, if bIMD is TRUE. */
2319 IMD_finalize(inputrec->bIMD, inputrec->imd);
2321 if (converged)
2323 step--; /* we never took that last step in this case */
2326 if (ems.fmax > inputrec->em_tol)
2328 if (MASTER(cr))
2330 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2331 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2333 converged = FALSE;
2336 /* If we printed energy and/or logfile last step (which was the last step)
2337 * we don't have to do it again, but otherwise print the final values.
2339 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2341 print_ebin_header(fplog, step, step);
2343 if (!do_ene || !do_log) /* Write final energy file entries */
2345 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2346 !do_log ? fplog : nullptr, step, step, eprNORMAL,
2347 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2350 /* Print some stuff... */
2351 if (MASTER(cr))
2353 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2356 /* IMPORTANT!
2357 * For accurate normal mode calculation it is imperative that we
2358 * store the last conformation into the full precision binary trajectory.
2360 * However, we should only do it if we did NOT already write this step
2361 * above (which we did if do_x or do_f was true).
2363 do_x = !do_per_step(step, inputrec->nstxout);
2364 do_f = !do_per_step(step, inputrec->nstfout);
2365 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2366 top_global, inputrec, step,
2367 &ems, state_global, observablesHistory);
2369 if (MASTER(cr))
2371 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2372 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2373 number_steps, &ems, sqrtNumAtoms);
2374 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2375 number_steps, &ems, sqrtNumAtoms);
2377 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2380 finish_em(cr, outf, walltime_accounting, wcycle);
2382 /* To print the actual number of steps we needed somewhere */
2383 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2385 return 0;
2386 } /* That's all folks */
2388 /*! \brief Do steepest descents minimization
2389 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
2390 int nfile, const t_filenm fnm[],
2391 const gmx_output_env_t *oenv,
2392 const MdrunOptions &mdrunOptions,
2393 gmx_vsite_t *vsite, gmx_constr_t constr,
2394 gmx::IMDOutputProvider *outputProvider,
2395 t_inputrec *inputrec,
2396 gmx_mtop_t *top_global, t_fcdata *fcd,
2397 t_state *state_global,
2398 gmx::MDAtoms *mdAtoms,
2399 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2400 gmx_edsam_t ed,
2401 t_forcerec *fr,
2402 const ReplicaExchangeParameters &replExParams,
2403 gmx_walltime_accounting_t walltime_accounting)
2405 double do_steep(FILE *fplog, t_commrec *cr, const gmx::MDLogger gmx_unused &mdlog,
2406 int nfile, const t_filenm fnm[],
2407 const gmx_output_env_t gmx_unused *oenv,
2408 const MdrunOptions &mdrunOptions,
2409 gmx_vsite_t *vsite, gmx_constr_t constr,
2410 gmx::IMDOutputProvider *outputProvider,
2411 t_inputrec *inputrec,
2412 gmx_mtop_t *top_global, t_fcdata *fcd,
2413 t_state *state_global,
2414 ObservablesHistory *observablesHistory,
2415 gmx::MDAtoms *mdAtoms,
2416 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2417 t_forcerec *fr,
2418 const ReplicaExchangeParameters gmx_unused &replExParams,
2419 gmx_membed_t gmx_unused *membed,
2420 gmx_walltime_accounting_t walltime_accounting)
2422 const char *SD = "Steepest Descents";
2423 gmx_localtop_t *top;
2424 gmx_enerdata_t *enerd;
2425 gmx_global_stat_t gstat;
2426 t_graph *graph;
2427 real stepsize;
2428 real ustep;
2429 gmx_mdoutf_t outf;
2430 t_mdebin *mdebin;
2431 gmx_bool bDone, bAbort, do_x, do_f;
2432 tensor vir, pres;
2433 rvec mu_tot;
2434 int nsteps;
2435 int count = 0;
2436 int steps_accepted = 0;
2437 auto mdatoms = mdAtoms->mdatoms();
2439 /* Create 2 states on the stack and extract pointers that we will swap */
2440 em_state_t s0 {}, s1 {};
2441 em_state_t *s_min = &s0;
2442 em_state_t *s_try = &s1;
2444 /* Init em and store the local state in s_try */
2445 init_em(fplog, SD, cr, outputProvider, inputrec, mdrunOptions,
2446 state_global, top_global, s_try, &top,
2447 nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
2448 vsite, constr, nullptr,
2449 nfile, fnm, &outf, &mdebin, wcycle);
2451 /* Print to log file */
2452 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2454 /* Set variables for stepsize (in nm). This is the largest
2455 * step that we are going to make in any direction.
2457 ustep = inputrec->em_stepsize;
2458 stepsize = 0;
2460 /* Max number of steps */
2461 nsteps = inputrec->nsteps;
2463 if (MASTER(cr))
2465 /* Print to the screen */
2466 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2468 if (fplog)
2470 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2473 /**** HERE STARTS THE LOOP ****
2474 * count is the counter for the number of steps
2475 * bDone will be TRUE when the minimization has converged
2476 * bAbort will be TRUE when nsteps steps have been performed or when
2477 * the stepsize becomes smaller than is reasonable for machine precision
2479 count = 0;
2480 bDone = FALSE;
2481 bAbort = FALSE;
2482 while (!bDone && !bAbort)
2484 bAbort = (nsteps >= 0) && (count == nsteps);
2486 /* set new coordinates, except for first step */
2487 bool validStep = true;
2488 if (count > 0)
2490 validStep =
2491 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2492 s_min, stepsize, &s_min->f, s_try,
2493 constr, top, nrnb, wcycle, count);
2496 if (validStep)
2498 evaluate_energy(fplog, cr,
2499 top_global, s_try, top,
2500 inputrec, nrnb, wcycle, gstat,
2501 vsite, constr, fcd, graph, mdAtoms, fr,
2502 mu_tot, enerd, vir, pres, count, count == 0);
2504 else
2506 // Signal constraint error during stepping with energy=inf
2507 s_try->epot = std::numeric_limits<real>::infinity();
2510 if (MASTER(cr))
2512 print_ebin_header(fplog, count, count);
2515 if (count == 0)
2517 s_min->epot = s_try->epot;
2520 /* Print it if necessary */
2521 if (MASTER(cr))
2523 if (mdrunOptions.verbose)
2525 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2526 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2527 ( (count == 0) || (s_try->epot < s_min->epot) ) ? '\n' : '\r');
2528 fflush(stderr);
2531 if ( (count == 0) || (s_try->epot < s_min->epot) )
2533 /* Store the new (lower) energies */
2534 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2535 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2536 s_try->s.box, nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2538 /* Prepare IMD energy record, if bIMD is TRUE. */
2539 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2541 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2542 do_per_step(steps_accepted, inputrec->nstdisreout),
2543 do_per_step(steps_accepted, inputrec->nstorireout),
2544 fplog, count, count, eprNORMAL,
2545 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2546 fflush(fplog);
2550 /* Now if the new energy is smaller than the previous...
2551 * or if this is the first step!
2552 * or if we did random steps!
2555 if ( (count == 0) || (s_try->epot < s_min->epot) )
2557 steps_accepted++;
2559 /* Test whether the convergence criterion is met... */
2560 bDone = (s_try->fmax < inputrec->em_tol);
2562 /* Copy the arrays for force, positions and energy */
2563 /* The 'Min' array always holds the coords and forces of the minimal
2564 sampled energy */
2565 swap_em_state(&s_min, &s_try);
2566 if (count > 0)
2568 ustep *= 1.2;
2571 /* Write to trn, if necessary */
2572 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2573 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2574 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
2575 top_global, inputrec, count,
2576 s_min, state_global, observablesHistory);
2578 else
2580 /* If energy is not smaller make the step smaller... */
2581 ustep *= 0.5;
2583 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2585 /* Reload the old state */
2586 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2587 s_min, top, mdAtoms, fr, vsite, constr,
2588 nrnb, wcycle);
2592 /* Determine new step */
2593 stepsize = ustep/s_min->fmax;
2595 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2596 #if GMX_DOUBLE
2597 if (count == nsteps || ustep < 1e-12)
2598 #else
2599 if (count == nsteps || ustep < 1e-6)
2600 #endif
2602 if (MASTER(cr))
2604 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != nullptr);
2605 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != nullptr);
2607 bAbort = TRUE;
2610 /* Send IMD energies and positions, if bIMD is TRUE. */
2611 if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box,
2612 MASTER(cr) ? as_rvec_array(state_global->x.data()) : nullptr,
2613 inputrec, 0, wcycle) &&
2614 MASTER(cr))
2616 IMD_send_positions(inputrec->imd);
2619 count++;
2620 } /* End of the loop */
2622 /* IMD cleanup, if bIMD is TRUE. */
2623 IMD_finalize(inputrec->bIMD, inputrec->imd);
2625 /* Print some data... */
2626 if (MASTER(cr))
2628 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2630 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2631 top_global, inputrec, count,
2632 s_min, state_global, observablesHistory);
2634 if (MASTER(cr))
2636 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2638 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2639 s_min, sqrtNumAtoms);
2640 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2641 s_min, sqrtNumAtoms);
2644 finish_em(cr, outf, walltime_accounting, wcycle);
2646 /* To print the actual number of steps we needed somewhere */
2647 inputrec->nsteps = count;
2649 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2651 return 0;
2652 } /* That's all folks */
2654 /*! \brief Do normal modes analysis
2655 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
2656 int nfile, const t_filenm fnm[],
2657 const gmx_output_env_t *oenv,
2658 const MdrunOptions &mdrunOptions,
2659 gmx_vsite_t *vsite, gmx_constr_t constr,
2660 gmx::IMDOutputProvider *outputProvider,
2661 t_inputrec *inputrec,
2662 gmx_mtop_t *top_global, t_fcdata *fcd,
2663 t_state *state_global,
2664 gmx::MDAtoms *mdAtoms,
2665 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2666 gmx_edsam_t ed,
2667 t_forcerec *fr,
2668 const ReplicaExchangeParameters &replExParams,
2669 gmx_walltime_accounting_t walltime_accounting)
2671 double do_nm(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
2672 int nfile, const t_filenm fnm[],
2673 const gmx_output_env_t gmx_unused *oenv,
2674 const MdrunOptions &mdrunOptions,
2675 gmx_vsite_t *vsite, gmx_constr_t constr,
2676 gmx::IMDOutputProvider *outputProvider,
2677 t_inputrec *inputrec,
2678 gmx_mtop_t *top_global, t_fcdata *fcd,
2679 t_state *state_global,
2680 ObservablesHistory gmx_unused *observablesHistory,
2681 gmx::MDAtoms *mdAtoms,
2682 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2683 t_forcerec *fr,
2684 const ReplicaExchangeParameters gmx_unused &replExParams,
2685 gmx_membed_t gmx_unused *membed,
2686 gmx_walltime_accounting_t walltime_accounting)
2688 const char *NM = "Normal Mode Analysis";
2689 gmx_mdoutf_t outf;
2690 int nnodes, node;
2691 gmx_localtop_t *top;
2692 gmx_enerdata_t *enerd;
2693 gmx_global_stat_t gstat;
2694 t_graph *graph;
2695 tensor vir, pres;
2696 rvec mu_tot;
2697 rvec *fneg, *dfdx;
2698 gmx_bool bSparse; /* use sparse matrix storage format */
2699 size_t sz;
2700 gmx_sparsematrix_t * sparse_matrix = nullptr;
2701 real * full_matrix = nullptr;
2703 /* added with respect to mdrun */
2704 int row, col;
2705 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2706 real x_min;
2707 bool bIsMaster = MASTER(cr);
2708 auto mdatoms = mdAtoms->mdatoms();
2710 if (constr != nullptr)
2712 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2715 gmx_shellfc_t *shellfc;
2717 em_state_t state_work {};
2719 /* Init em and store the local state in state_minimum */
2720 init_em(fplog, NM, cr, outputProvider, inputrec, mdrunOptions,
2721 state_global, top_global, &state_work, &top,
2722 nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
2723 vsite, constr, &shellfc,
2724 nfile, fnm, &outf, nullptr, wcycle);
2726 std::vector<size_t> atom_index = get_atom_index(top_global);
2727 snew(fneg, atom_index.size());
2728 snew(dfdx, atom_index.size());
2730 #if !GMX_DOUBLE
2731 if (bIsMaster)
2733 fprintf(stderr,
2734 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2735 " which MIGHT not be accurate enough for normal mode analysis.\n"
2736 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
2737 " are fairly modest even if you recompile in double precision.\n\n");
2739 #endif
2741 /* Check if we can/should use sparse storage format.
2743 * Sparse format is only useful when the Hessian itself is sparse, which it
2744 * will be when we use a cutoff.
2745 * For small systems (n<1000) it is easier to always use full matrix format, though.
2747 if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
2749 GMX_LOG(mdlog.warning).appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2750 bSparse = FALSE;
2752 else if (atom_index.size() < 1000)
2754 GMX_LOG(mdlog.warning).appendTextFormatted("Small system size (N=%d), using full Hessian format.",
2755 atom_index.size());
2756 bSparse = FALSE;
2758 else
2760 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2761 bSparse = TRUE;
2764 /* Number of dimensions, based on real atoms, that is not vsites or shell */
2765 sz = DIM*atom_index.size();
2767 fprintf(stderr, "Allocating Hessian memory...\n\n");
2769 if (bSparse)
2771 sparse_matrix = gmx_sparsematrix_init(sz);
2772 sparse_matrix->compressed_symmetric = TRUE;
2774 else
2776 snew(full_matrix, sz*sz);
2779 init_nrnb(nrnb);
2781 where();
2783 /* Write start time and temperature */
2784 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2786 /* fudge nr of steps to nr of atoms */
2787 inputrec->nsteps = atom_index.size()*2;
2789 if (bIsMaster)
2791 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2792 *(top_global->name), (int)inputrec->nsteps);
2795 nnodes = cr->nnodes;
2797 /* Make evaluate_energy do a single node force calculation */
2798 cr->nnodes = 1;
2799 evaluate_energy(fplog, cr,
2800 top_global, &state_work, top,
2801 inputrec, nrnb, wcycle, gstat,
2802 vsite, constr, fcd, graph, mdAtoms, fr,
2803 mu_tot, enerd, vir, pres, -1, TRUE);
2804 cr->nnodes = nnodes;
2806 /* if forces are not small, warn user */
2807 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
2809 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
2810 if (state_work.fmax > 1.0e-3)
2812 GMX_LOG(mdlog.warning).appendText(
2813 "The force is probably not small enough to "
2814 "ensure that you are at a minimum.\n"
2815 "Be aware that negative eigenvalues may occur\n"
2816 "when the resulting matrix is diagonalized.");
2819 /***********************************************************
2821 * Loop over all pairs in matrix
2823 * do_force called twice. Once with positive and
2824 * once with negative displacement
2826 ************************************************************/
2828 /* Steps are divided one by one over the nodes */
2829 bool bNS = true;
2830 for (unsigned int aid = cr->nodeid; aid < atom_index.size(); aid += nnodes)
2832 size_t atom = atom_index[aid];
2833 for (size_t d = 0; d < DIM; d++)
2835 gmx_bool bBornRadii = FALSE;
2836 gmx_int64_t step = 0;
2837 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2838 double t = 0;
2840 x_min = state_work.s.x[atom][d];
2842 for (unsigned int dx = 0; (dx < 2); dx++)
2844 if (dx == 0)
2846 state_work.s.x[atom][d] = x_min - der_range;
2848 else
2850 state_work.s.x[atom][d] = x_min + der_range;
2853 /* Make evaluate_energy do a single node force calculation */
2854 cr->nnodes = 1;
2855 if (shellfc)
2857 /* Now is the time to relax the shells */
2858 (void) relax_shell_flexcon(fplog, cr, mdrunOptions.verbose, step,
2859 inputrec, bNS, force_flags,
2860 top,
2861 constr, enerd, fcd,
2862 &state_work.s, &state_work.f, vir, mdatoms,
2863 nrnb, wcycle, graph, &top_global->groups,
2864 shellfc, fr, bBornRadii, t, mu_tot,
2865 vsite,
2866 DdOpenBalanceRegionBeforeForceComputation::no,
2867 DdCloseBalanceRegionAfterForceComputation::no);
2868 bNS = false;
2869 step++;
2871 else
2873 evaluate_energy(fplog, cr,
2874 top_global, &state_work, top,
2875 inputrec, nrnb, wcycle, gstat,
2876 vsite, constr, fcd, graph, mdAtoms, fr,
2877 mu_tot, enerd, vir, pres, atom*2+dx, FALSE);
2880 cr->nnodes = nnodes;
2882 if (dx == 0)
2884 for (size_t i = 0; i < atom_index.size(); i++)
2886 copy_rvec(state_work.f[atom_index[i]], fneg[i]);
2891 /* x is restored to original */
2892 state_work.s.x[atom][d] = x_min;
2894 for (size_t j = 0; j < atom_index.size(); j++)
2896 for (size_t k = 0; (k < DIM); k++)
2898 dfdx[j][k] =
2899 -(state_work.f[atom_index[j]][k] - fneg[j][k])/(2*der_range);
2903 if (!bIsMaster)
2905 #if GMX_MPI
2906 #define mpi_type GMX_MPI_REAL
2907 MPI_Send(dfdx[0], atom_index.size()*DIM, mpi_type, MASTER(cr),
2908 cr->nodeid, cr->mpi_comm_mygroup);
2909 #endif
2911 else
2913 for (node = 0; (node < nnodes && atom+node < atom_index.size()); node++)
2915 if (node > 0)
2917 #if GMX_MPI
2918 MPI_Status stat;
2919 MPI_Recv(dfdx[0], atom_index.size()*DIM, mpi_type, node, node,
2920 cr->mpi_comm_mygroup, &stat);
2921 #undef mpi_type
2922 #endif
2925 row = (atom + node)*DIM + d;
2927 for (size_t j = 0; j < atom_index.size(); j++)
2929 for (size_t k = 0; k < DIM; k++)
2931 col = j*DIM + k;
2933 if (bSparse)
2935 if (col >= row && dfdx[j][k] != 0.0)
2937 gmx_sparsematrix_increment_value(sparse_matrix,
2938 row, col, dfdx[j][k]);
2941 else
2943 full_matrix[row*sz+col] = dfdx[j][k];
2950 if (mdrunOptions.verbose && fplog)
2952 fflush(fplog);
2955 /* write progress */
2956 if (bIsMaster && mdrunOptions.verbose)
2958 fprintf(stderr, "\rFinished step %d out of %d",
2959 static_cast<int>(std::min(atom+nnodes, atom_index.size())),
2960 static_cast<int>(atom_index.size()));
2961 fflush(stderr);
2965 if (bIsMaster)
2967 fprintf(stderr, "\n\nWriting Hessian...\n");
2968 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2971 finish_em(cr, outf, walltime_accounting, wcycle);
2973 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size()*2);
2975 return 0;
2978 } // namespace gmx