Introduce ObservablesHistory container
[gromacs.git] / src / gromacs / mdlib / minimize.cpp
blob55937087f7e52672b7b317caa77d68db946d6757
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 t_state *state_global, gmx_mtop_t *top_global,
328 em_state_t *ems, gmx_localtop_t **top,
329 t_nrnb *nrnb, rvec mu_tot,
330 t_forcerec *fr, gmx_enerdata_t **enerd,
331 t_graph **graph, t_mdatoms *mdatoms, gmx_global_stat_t *gstat,
332 gmx_vsite_t *vsite, gmx_constr_t constr, gmx_shellfc_t **shellfc,
333 int nfile, const t_filenm fnm[],
334 gmx_mdoutf_t *outf, t_mdebin **mdebin,
335 int imdport, unsigned long gmx_unused Flags,
336 gmx_wallcycle_t wcycle)
338 real dvdl_constr;
340 if (fplog)
342 fprintf(fplog, "Initiating %s\n", title);
345 state_global->ngtc = 0;
347 /* Initialize lambda variables */
348 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, nullptr);
350 init_nrnb(nrnb);
352 /* Interactive molecular dynamics */
353 init_IMD(ir, cr, top_global, fplog, 1, as_rvec_array(state_global->x.data()),
354 nfile, fnm, nullptr, imdport, Flags);
356 if (ir->eI == eiNM)
358 GMX_ASSERT(shellfc != NULL, "With NM we always support shells");
360 *shellfc = init_shell_flexcon(stdout,
361 top_global,
362 n_flexible_constraints(constr),
363 ir->nstcalcenergy,
364 DOMAINDECOMP(cr));
366 else
368 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI), "This else currently only handles energy minimizers, consider if your algorithm needs shell/flexible-constraint support");
370 /* With energy minimization, shells and flexible constraints are
371 * automatically minimized when treated like normal DOFS.
373 if (shellfc != nullptr)
375 *shellfc = nullptr;
379 if (DOMAINDECOMP(cr))
381 *top = dd_init_local_top(top_global);
383 dd_init_local_state(cr->dd, state_global, &ems->s);
385 /* Distribute the charge groups over the nodes from the master node */
386 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
387 state_global, top_global, ir,
388 &ems->s, &ems->f, mdatoms, *top,
389 fr, vsite, constr,
390 nrnb, nullptr, FALSE);
391 dd_store_state(cr->dd, &ems->s);
393 *graph = nullptr;
395 else
397 state_change_natoms(state_global, state_global->natoms);
398 /* Just copy the state */
399 ems->s = *state_global;
400 state_change_natoms(&ems->s, ems->s.natoms);
401 /* We need to allocate one element extra, since we might use
402 * (unaligned) 4-wide SIMD loads to access rvec entries.
404 ems->f.resize(ems->s.natoms + 1);
406 snew(*top, 1);
407 mdAlgorithmsSetupAtomData(cr, ir, top_global, *top, fr,
408 graph, mdatoms,
409 vsite, shellfc ? *shellfc : nullptr);
411 if (vsite)
413 set_vsite_top(vsite, *top, mdatoms, cr);
417 update_mdatoms(mdatoms, state_global->lambda[efptMASS]);
419 if (constr)
421 if (ir->eConstrAlg == econtSHAKE &&
422 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
424 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
425 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
428 if (!DOMAINDECOMP(cr))
430 set_constraints(constr, *top, ir, mdatoms, cr);
433 if (!ir->bContinuation)
435 /* Constrain the starting coordinates */
436 dvdl_constr = 0;
437 constrain(PAR(cr) ? nullptr : fplog, TRUE, TRUE, constr, &(*top)->idef,
438 ir, cr, -1, 0, 1.0, mdatoms,
439 as_rvec_array(ems->s.x.data()),
440 as_rvec_array(ems->s.x.data()),
441 nullptr,
442 fr->bMolPBC, ems->s.box,
443 ems->s.lambda[efptFEP], &dvdl_constr,
444 nullptr, nullptr, nrnb, econqCoord);
448 if (PAR(cr))
450 *gstat = global_stat_init(ir);
452 else
454 *gstat = nullptr;
457 *outf = init_mdoutf(fplog, nfile, fnm, 0, cr, outputProvider, ir, top_global, nullptr, wcycle);
459 snew(*enerd, 1);
460 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
461 *enerd);
463 if (mdebin != nullptr)
465 /* Init bin for energy stuff */
466 *mdebin = init_mdebin(mdoutf_get_fp_ene(*outf), top_global, ir, nullptr);
469 clear_rvec(mu_tot);
470 calc_shifts(ems->s.box, fr->shift_vec);
473 //! Finalize the minimization
474 static void finish_em(t_commrec *cr, gmx_mdoutf_t outf,
475 gmx_walltime_accounting_t walltime_accounting,
476 gmx_wallcycle_t wcycle)
478 if (!(cr->duty & DUTY_PME))
480 /* Tell the PME only node to finish */
481 gmx_pme_send_finish(cr);
484 done_mdoutf(outf);
486 em_time_end(walltime_accounting, wcycle);
489 //! Swap two different EM states during minimization
490 static void swap_em_state(em_state_t **ems1, em_state_t **ems2)
492 em_state_t *tmp;
494 tmp = *ems1;
495 *ems1 = *ems2;
496 *ems2 = tmp;
499 //! Save the EM trajectory
500 static void write_em_traj(FILE *fplog, t_commrec *cr,
501 gmx_mdoutf_t outf,
502 gmx_bool bX, gmx_bool bF, const char *confout,
503 gmx_mtop_t *top_global,
504 t_inputrec *ir, gmx_int64_t step,
505 em_state_t *state,
506 t_state *state_global,
507 ObservablesHistory *observablesHistory)
509 int mdof_flags = 0;
511 if (bX)
513 mdof_flags |= MDOF_X;
515 if (bF)
517 mdof_flags |= MDOF_F;
520 /* If we want IMD output, set appropriate MDOF flag */
521 if (ir->bIMD)
523 mdof_flags |= MDOF_IMD;
526 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
527 top_global, step, (double)step,
528 &state->s, state_global, observablesHistory,
529 &state->f);
531 if (confout != nullptr && MASTER(cr))
533 GMX_RELEASE_ASSERT(bX, "The code below assumes that (with domain decomposition), x is collected to state_global in the call above.");
534 /* With domain decomposition the call above collected the state->s.x
535 * into state_global->x. Without DD we copy the local state pointer.
537 if (!DOMAINDECOMP(cr))
539 state_global = &state->s;
542 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
544 /* Make molecules whole only for confout writing */
545 do_pbc_mtop(fplog, ir->ePBC, state->s.box, top_global,
546 as_rvec_array(state_global->x.data()));
549 write_sto_conf_mtop(confout,
550 *top_global->name, top_global,
551 as_rvec_array(state_global->x.data()), nullptr, ir->ePBC, state->s.box);
555 //! \brief Do one minimization step
557 // \returns true when the step succeeded, false when a constraint error occurred
558 static bool do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
559 gmx_bool bMolPBC,
560 em_state_t *ems1, real a, const PaddedRVecVector *force,
561 em_state_t *ems2,
562 gmx_constr_t constr, gmx_localtop_t *top,
563 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
564 gmx_int64_t count)
567 t_state *s1, *s2;
568 int start, end;
569 real dvdl_constr;
570 int nthreads gmx_unused;
572 bool validStep = true;
574 s1 = &ems1->s;
575 s2 = &ems2->s;
577 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
579 gmx_incons("state mismatch in do_em_step");
582 s2->flags = s1->flags;
584 if (s2->natoms != s1->natoms)
586 state_change_natoms(s2, s1->natoms);
587 /* We need to allocate one element extra, since we might use
588 * (unaligned) 4-wide SIMD loads to access rvec entries.
590 ems2->f.resize(s2->natoms + 1);
592 if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
594 s2->cg_gl.resize(s1->cg_gl.size());
597 copy_mat(s1->box, s2->box);
598 /* Copy free energy state */
599 s2->lambda = s1->lambda;
600 copy_mat(s1->box, s2->box);
602 start = 0;
603 end = md->homenr;
605 // cppcheck-suppress unreadVariable
606 nthreads = gmx_omp_nthreads_get(emntUpdate);
607 #pragma omp parallel num_threads(nthreads)
609 const rvec *x1 = as_rvec_array(s1->x.data());
610 rvec *x2 = as_rvec_array(s2->x.data());
611 const rvec *f = as_rvec_array(force->data());
613 int gf = 0;
614 #pragma omp for schedule(static) nowait
615 for (int i = start; i < end; i++)
619 if (md->cFREEZE)
621 gf = md->cFREEZE[i];
623 for (int m = 0; m < DIM; m++)
625 if (ir->opts.nFreeze[gf][m])
627 x2[i][m] = x1[i][m];
629 else
631 x2[i][m] = x1[i][m] + a*f[i][m];
635 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
638 if (s2->flags & (1<<estCGP))
640 /* Copy the CG p vector */
641 const rvec *p1 = as_rvec_array(s1->cg_p.data());
642 rvec *p2 = as_rvec_array(s2->cg_p.data());
643 #pragma omp for schedule(static) nowait
644 for (int i = start; i < end; i++)
646 // Trivial OpenMP block that does not throw
647 copy_rvec(p1[i], p2[i]);
651 if (DOMAINDECOMP(cr))
653 s2->ddp_count = s1->ddp_count;
655 /* OpenMP does not supported unsigned loop variables */
656 #pragma omp for schedule(static) nowait
657 for (int i = 0; i < static_cast<int>(s2->cg_gl.size()); i++)
659 s2->cg_gl[i] = s1->cg_gl[i];
661 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
665 if (constr)
667 wallcycle_start(wcycle, ewcCONSTR);
668 dvdl_constr = 0;
669 validStep =
670 constrain(nullptr, TRUE, TRUE, constr, &top->idef,
671 ir, cr, count, 0, 1.0, md,
672 as_rvec_array(s1->x.data()), as_rvec_array(s2->x.data()),
673 nullptr, bMolPBC, s2->box,
674 s2->lambda[efptBONDED], &dvdl_constr,
675 nullptr, nullptr, nrnb, econqCoord);
676 wallcycle_stop(wcycle, ewcCONSTR);
678 // We should move this check to the different minimizers
679 if (!validStep && ir->eI != eiSteep)
681 gmx_fatal(FARGS, "The coordinates could not be constrained. Minimizer '%s' can not handle constraint failures, use minimizer '%s' before using '%s'.",
682 EI(ir->eI), EI(eiSteep), EI(ir->eI));
686 return validStep;
689 //! Prepare EM for using domain decomposition parallellization
690 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
691 gmx_mtop_t *top_global, t_inputrec *ir,
692 em_state_t *ems, gmx_localtop_t *top,
693 t_mdatoms *mdatoms, t_forcerec *fr,
694 gmx_vsite_t *vsite, gmx_constr_t constr,
695 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
697 /* Repartition the domain decomposition */
698 dd_partition_system(fplog, step, cr, FALSE, 1,
699 nullptr, top_global, ir,
700 &ems->s, &ems->f,
701 mdatoms, top, fr, vsite, constr,
702 nrnb, wcycle, FALSE);
703 dd_store_state(cr->dd, &ems->s);
706 //! De one energy evaluation
707 static void evaluate_energy(FILE *fplog, t_commrec *cr,
708 gmx_mtop_t *top_global,
709 em_state_t *ems, gmx_localtop_t *top,
710 t_inputrec *inputrec,
711 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
712 gmx_global_stat_t gstat,
713 gmx_vsite_t *vsite, gmx_constr_t constr,
714 t_fcdata *fcd,
715 t_graph *graph, t_mdatoms *mdatoms,
716 t_forcerec *fr, rvec mu_tot,
717 gmx_enerdata_t *enerd, tensor vir, tensor pres,
718 gmx_int64_t count, gmx_bool bFirst)
720 real t;
721 gmx_bool bNS;
722 tensor force_vir, shake_vir, ekin;
723 real dvdl_constr, prescorr, enercorr, dvdlcorr;
724 real terminate = 0;
726 /* Set the time to the initial time, the time does not change during EM */
727 t = inputrec->init_t;
729 if (bFirst ||
730 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
732 /* This is the first state or an old state used before the last ns */
733 bNS = TRUE;
735 else
737 bNS = FALSE;
738 if (inputrec->nstlist > 0)
740 bNS = TRUE;
744 if (vsite)
746 construct_vsites(vsite, as_rvec_array(ems->s.x.data()), 1, nullptr,
747 top->idef.iparams, top->idef.il,
748 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
751 if (DOMAINDECOMP(cr) && bNS)
753 /* Repartition the domain decomposition */
754 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
755 ems, top, mdatoms, fr, vsite, constr,
756 nrnb, wcycle);
759 /* Calc force & energy on new trial position */
760 /* do_force always puts the charge groups in the box and shifts again
761 * We do not unshift, so molecules are always whole in congrad.c
763 do_force(fplog, cr, inputrec,
764 count, nrnb, wcycle, top, &top_global->groups,
765 ems->s.box, &ems->s.x, &ems->s.hist,
766 &ems->f, force_vir, mdatoms, enerd, fcd,
767 ems->s.lambda, graph, fr, vsite, mu_tot, t, nullptr, TRUE,
768 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
769 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
770 (bNS ? GMX_FORCE_NS : 0));
772 /* Clear the unused shake virial and pressure */
773 clear_mat(shake_vir);
774 clear_mat(pres);
776 /* Communicate stuff when parallel */
777 if (PAR(cr) && inputrec->eI != eiNM)
779 wallcycle_start(wcycle, ewcMoveE);
781 global_stat(gstat, cr, enerd, force_vir, shake_vir, mu_tot,
782 inputrec, nullptr, nullptr, nullptr, 1, &terminate,
783 nullptr, FALSE,
784 CGLO_ENERGY |
785 CGLO_PRESSURE |
786 CGLO_CONSTRAINT);
788 wallcycle_stop(wcycle, ewcMoveE);
791 /* Calculate long range corrections to pressure and energy */
792 calc_dispcorr(inputrec, fr, ems->s.box, ems->s.lambda[efptVDW],
793 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
794 enerd->term[F_DISPCORR] = enercorr;
795 enerd->term[F_EPOT] += enercorr;
796 enerd->term[F_PRES] += prescorr;
797 enerd->term[F_DVDL] += dvdlcorr;
799 ems->epot = enerd->term[F_EPOT];
801 if (constr)
803 /* Project out the constraint components of the force */
804 wallcycle_start(wcycle, ewcCONSTR);
805 dvdl_constr = 0;
806 rvec *f_rvec = as_rvec_array(ems->f.data());
807 constrain(nullptr, FALSE, FALSE, constr, &top->idef,
808 inputrec, cr, count, 0, 1.0, mdatoms,
809 as_rvec_array(ems->s.x.data()), f_rvec, f_rvec,
810 fr->bMolPBC, ems->s.box,
811 ems->s.lambda[efptBONDED], &dvdl_constr,
812 nullptr, &shake_vir, nrnb, econqForceDispl);
813 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
814 m_add(force_vir, shake_vir, vir);
815 wallcycle_stop(wcycle, ewcCONSTR);
817 else
819 copy_mat(force_vir, vir);
822 clear_mat(ekin);
823 enerd->term[F_PRES] =
824 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
826 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
828 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
830 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
834 //! Parallel utility summing energies and forces
835 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
836 gmx_mtop_t *top_global,
837 em_state_t *s_min, em_state_t *s_b)
839 t_block *cgs_gl;
840 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
841 double partsum;
842 unsigned char *grpnrFREEZE;
844 if (debug)
846 fprintf(debug, "Doing reorder_partsum\n");
849 const rvec *fm = as_rvec_array(s_min->f.data());
850 const rvec *fb = as_rvec_array(s_b->f.data());
852 cgs_gl = dd_charge_groups_global(cr->dd);
853 index = cgs_gl->index;
855 /* Collect fm in a global vector fmg.
856 * This conflicts with the spirit of domain decomposition,
857 * but to fully optimize this a much more complicated algorithm is required.
859 rvec *fmg;
860 snew(fmg, top_global->natoms);
862 ncg = s_min->s.cg_gl.size();
863 cg_gl = s_min->s.cg_gl.data();
864 i = 0;
865 for (c = 0; c < ncg; c++)
867 cg = cg_gl[c];
868 a0 = index[cg];
869 a1 = index[cg+1];
870 for (a = a0; a < a1; a++)
872 copy_rvec(fm[i], fmg[a]);
873 i++;
876 gmx_sum(top_global->natoms*3, fmg[0], cr);
878 /* Now we will determine the part of the sum for the cgs in state s_b */
879 ncg = s_b->s.cg_gl.size();
880 cg_gl = s_b->s.cg_gl.data();
881 partsum = 0;
882 i = 0;
883 gf = 0;
884 grpnrFREEZE = top_global->groups.grpnr[egcFREEZE];
885 for (c = 0; c < ncg; c++)
887 cg = cg_gl[c];
888 a0 = index[cg];
889 a1 = index[cg+1];
890 for (a = a0; a < a1; a++)
892 if (mdatoms->cFREEZE && grpnrFREEZE)
894 gf = grpnrFREEZE[i];
896 for (m = 0; m < DIM; m++)
898 if (!opts->nFreeze[gf][m])
900 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
903 i++;
907 sfree(fmg);
909 return partsum;
912 //! Print some stuff, like beta, whatever that means.
913 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
914 gmx_mtop_t *top_global,
915 em_state_t *s_min, em_state_t *s_b)
917 double sum;
919 /* This is just the classical Polak-Ribiere calculation of beta;
920 * it looks a bit complicated since we take freeze groups into account,
921 * and might have to sum it in parallel runs.
924 if (!DOMAINDECOMP(cr) ||
925 (s_min->s.ddp_count == cr->dd->ddp_count &&
926 s_b->s.ddp_count == cr->dd->ddp_count))
928 const rvec *fm = as_rvec_array(s_min->f.data());
929 const rvec *fb = as_rvec_array(s_b->f.data());
930 sum = 0;
931 int gf = 0;
932 /* This part of code can be incorrect with DD,
933 * since the atom ordering in s_b and s_min might differ.
935 for (int i = 0; i < mdatoms->homenr; i++)
937 if (mdatoms->cFREEZE)
939 gf = mdatoms->cFREEZE[i];
941 for (int m = 0; m < DIM; m++)
943 if (!opts->nFreeze[gf][m])
945 sum += (fb[i][m] - fm[i][m])*fb[i][m];
950 else
952 /* We need to reorder cgs while summing */
953 sum = reorder_partsum(cr, opts, mdatoms, top_global, s_min, s_b);
955 if (PAR(cr))
957 gmx_sumd(1, &sum, cr);
960 return sum/gmx::square(s_min->fnorm);
963 namespace gmx
966 /*! \brief Do conjugate gradients minimization
967 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
968 int nfile, const t_filenm fnm[],
969 const gmx_output_env_t *oenv, gmx_bool bVerbose,
970 int nstglobalcomm,
971 gmx_vsite_t *vsite, gmx_constr_t constr,
972 int stepout,
973 gmx::IMDOutputProvider *outputProvider,
974 t_inputrec *inputrec,
975 gmx_mtop_t *top_global, t_fcdata *fcd,
976 t_state *state_global,
977 t_mdatoms *mdatoms,
978 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
979 gmx_edsam_t ed,
980 t_forcerec *fr,
981 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
982 gmx_membed_t gmx_unused *membed,
983 real cpt_period, real max_hours,
984 int imdport,
985 unsigned long Flags,
986 gmx_walltime_accounting_t walltime_accounting)
988 double do_cg(FILE *fplog, t_commrec *cr, const gmx::MDLogger gmx_unused &mdlog,
989 int nfile, const t_filenm fnm[],
990 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
991 int gmx_unused nstglobalcomm,
992 gmx_vsite_t *vsite, gmx_constr_t constr,
993 int gmx_unused stepout,
994 gmx::IMDOutputProvider *outputProvider,
995 t_inputrec *inputrec,
996 gmx_mtop_t *top_global, t_fcdata *fcd,
997 t_state *state_global,
998 ObservablesHistory *observablesHistory,
999 t_mdatoms *mdatoms,
1000 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1001 gmx_edsam_t gmx_unused ed,
1002 t_forcerec *fr,
1003 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1004 gmx_membed_t gmx_unused *membed,
1005 real gmx_unused cpt_period, real gmx_unused max_hours,
1006 int imdport,
1007 unsigned long gmx_unused Flags,
1008 gmx_walltime_accounting_t walltime_accounting)
1010 const char *CG = "Polak-Ribiere Conjugate Gradients";
1012 gmx_localtop_t *top;
1013 gmx_enerdata_t *enerd;
1014 gmx_global_stat_t gstat;
1015 t_graph *graph;
1016 double tmp, minstep;
1017 real stepsize;
1018 real a, b, c, beta = 0.0;
1019 real epot_repl = 0;
1020 real pnorm;
1021 t_mdebin *mdebin;
1022 gmx_bool converged, foundlower;
1023 rvec mu_tot;
1024 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1025 tensor vir, pres;
1026 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1027 gmx_mdoutf_t outf;
1028 int m, step, nminstep;
1030 step = 0;
1032 // Ensure the extra per-atom state array gets allocated
1033 state_global->flags |= (1<<estCGP);
1035 /* Create 4 states on the stack and extract pointers that we will swap */
1036 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1037 em_state_t *s_min = &s0;
1038 em_state_t *s_a = &s1;
1039 em_state_t *s_b = &s2;
1040 em_state_t *s_c = &s3;
1042 /* Init em and store the local state in s_min */
1043 init_em(fplog, CG, cr, outputProvider, inputrec,
1044 state_global, top_global, s_min, &top,
1045 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat,
1046 vsite, constr, nullptr,
1047 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
1049 /* Print to log file */
1050 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1052 /* Max number of steps */
1053 number_steps = inputrec->nsteps;
1055 if (MASTER(cr))
1057 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1059 if (fplog)
1061 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1064 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1065 /* do_force always puts the charge groups in the box and shifts again
1066 * We do not unshift, so molecules are always whole in congrad.c
1068 evaluate_energy(fplog, cr,
1069 top_global, s_min, top,
1070 inputrec, nrnb, wcycle, gstat,
1071 vsite, constr, fcd, graph, mdatoms, fr,
1072 mu_tot, enerd, vir, pres, -1, TRUE);
1073 where();
1075 if (MASTER(cr))
1077 /* Copy stuff to the energy bin for easy printing etc. */
1078 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1079 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1080 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1082 print_ebin_header(fplog, step, step);
1083 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1084 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1086 where();
1088 /* Estimate/guess the initial stepsize */
1089 stepsize = inputrec->em_stepsize/s_min->fnorm;
1091 if (MASTER(cr))
1093 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1094 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1095 s_min->fmax, s_min->a_fmax+1);
1096 fprintf(stderr, " F-Norm = %12.5e\n",
1097 s_min->fnorm/sqrtNumAtoms);
1098 fprintf(stderr, "\n");
1099 /* and copy to the log file too... */
1100 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1101 s_min->fmax, s_min->a_fmax+1);
1102 fprintf(fplog, " F-Norm = %12.5e\n",
1103 s_min->fnorm/sqrtNumAtoms);
1104 fprintf(fplog, "\n");
1106 /* Start the loop over CG steps.
1107 * Each successful step is counted, and we continue until
1108 * we either converge or reach the max number of steps.
1110 converged = FALSE;
1111 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1114 /* start taking steps in a new direction
1115 * First time we enter the routine, beta=0, and the direction is
1116 * simply the negative gradient.
1119 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1120 rvec *pm = as_rvec_array(s_min->s.cg_p.data());
1121 const rvec *sfm = as_rvec_array(s_min->f.data());
1122 double gpa = 0;
1123 int gf = 0;
1124 for (int i = 0; i < mdatoms->homenr; i++)
1126 if (mdatoms->cFREEZE)
1128 gf = mdatoms->cFREEZE[i];
1130 for (m = 0; m < DIM; m++)
1132 if (!inputrec->opts.nFreeze[gf][m])
1134 pm[i][m] = sfm[i][m] + beta*pm[i][m];
1135 gpa -= pm[i][m]*sfm[i][m];
1136 /* f is negative gradient, thus the sign */
1138 else
1140 pm[i][m] = 0;
1145 /* Sum the gradient along the line across CPUs */
1146 if (PAR(cr))
1148 gmx_sumd(1, &gpa, cr);
1151 /* Calculate the norm of the search vector */
1152 get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1154 /* Just in case stepsize reaches zero due to numerical precision... */
1155 if (stepsize <= 0)
1157 stepsize = inputrec->em_stepsize/pnorm;
1161 * Double check the value of the derivative in the search direction.
1162 * If it is positive it must be due to the old information in the
1163 * CG formula, so just remove that and start over with beta=0.
1164 * This corresponds to a steepest descent step.
1166 if (gpa > 0)
1168 beta = 0;
1169 step--; /* Don't count this step since we are restarting */
1170 continue; /* Go back to the beginning of the big for-loop */
1173 /* Calculate minimum allowed stepsize, before the average (norm)
1174 * relative change in coordinate is smaller than precision
1176 minstep = 0;
1177 for (int i = 0; i < mdatoms->homenr; i++)
1179 for (m = 0; m < DIM; m++)
1181 tmp = fabs(s_min->s.x[i][m]);
1182 if (tmp < 1.0)
1184 tmp = 1.0;
1186 tmp = pm[i][m]/tmp;
1187 minstep += tmp*tmp;
1190 /* Add up from all CPUs */
1191 if (PAR(cr))
1193 gmx_sumd(1, &minstep, cr);
1196 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1198 if (stepsize < minstep)
1200 converged = TRUE;
1201 break;
1204 /* Write coordinates if necessary */
1205 do_x = do_per_step(step, inputrec->nstxout);
1206 do_f = do_per_step(step, inputrec->nstfout);
1208 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
1209 top_global, inputrec, step,
1210 s_min, state_global, observablesHistory);
1212 /* Take a step downhill.
1213 * In theory, we should minimize the function along this direction.
1214 * That is quite possible, but it turns out to take 5-10 function evaluations
1215 * for each line. However, we dont really need to find the exact minimum -
1216 * it is much better to start a new CG step in a modified direction as soon
1217 * as we are close to it. This will save a lot of energy evaluations.
1219 * In practice, we just try to take a single step.
1220 * If it worked (i.e. lowered the energy), we increase the stepsize but
1221 * the continue straight to the next CG step without trying to find any minimum.
1222 * If it didn't work (higher energy), there must be a minimum somewhere between
1223 * the old position and the new one.
1225 * Due to the finite numerical accuracy, it turns out that it is a good idea
1226 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1227 * This leads to lower final energies in the tests I've done. / Erik
1229 s_a->epot = s_min->epot;
1230 a = 0.0;
1231 c = a + stepsize; /* reference position along line is zero */
1233 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1235 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1236 s_min, top, mdatoms, fr, vsite, constr,
1237 nrnb, wcycle);
1240 /* Take a trial step (new coords in s_c) */
1241 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, &s_min->s.cg_p, s_c,
1242 constr, top, nrnb, wcycle, -1);
1244 neval++;
1245 /* Calculate energy for the trial step */
1246 evaluate_energy(fplog, cr,
1247 top_global, s_c, top,
1248 inputrec, nrnb, wcycle, gstat,
1249 vsite, constr, fcd, graph, mdatoms, fr,
1250 mu_tot, enerd, vir, pres, -1, FALSE);
1252 /* Calc derivative along line */
1253 const rvec *pc = as_rvec_array(s_c->s.cg_p.data());
1254 const rvec *sfc = as_rvec_array(s_c->f.data());
1255 double gpc = 0;
1256 for (int i = 0; i < mdatoms->homenr; i++)
1258 for (m = 0; m < DIM; m++)
1260 gpc -= pc[i][m]*sfc[i][m]; /* f is negative gradient, thus the sign */
1263 /* Sum the gradient along the line across CPUs */
1264 if (PAR(cr))
1266 gmx_sumd(1, &gpc, cr);
1269 /* This is the max amount of increase in energy we tolerate */
1270 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1272 /* Accept the step if the energy is lower, or if it is not significantly higher
1273 * and the line derivative is still negative.
1275 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1277 foundlower = TRUE;
1278 /* Great, we found a better energy. Increase step for next iteration
1279 * if we are still going down, decrease it otherwise
1281 if (gpc < 0)
1283 stepsize *= 1.618034; /* The golden section */
1285 else
1287 stepsize *= 0.618034; /* 1/golden section */
1290 else
1292 /* New energy is the same or higher. We will have to do some work
1293 * to find a smaller value in the interval. Take smaller step next time!
1295 foundlower = FALSE;
1296 stepsize *= 0.618034;
1302 /* OK, if we didn't find a lower value we will have to locate one now - there must
1303 * be one in the interval [a=0,c].
1304 * The same thing is valid here, though: Don't spend dozens of iterations to find
1305 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1306 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1308 * I also have a safeguard for potentially really pathological functions so we never
1309 * take more than 20 steps before we give up ...
1311 * If we already found a lower value we just skip this step and continue to the update.
1313 double gpb;
1314 if (!foundlower)
1316 nminstep = 0;
1320 /* Select a new trial point.
1321 * If the derivatives at points a & c have different sign we interpolate to zero,
1322 * otherwise just do a bisection.
1324 if (gpa < 0 && gpc > 0)
1326 b = a + gpa*(a-c)/(gpc-gpa);
1328 else
1330 b = 0.5*(a+c);
1333 /* safeguard if interpolation close to machine accuracy causes errors:
1334 * never go outside the interval
1336 if (b <= a || b >= c)
1338 b = 0.5*(a+c);
1341 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1343 /* Reload the old state */
1344 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1345 s_min, top, mdatoms, fr, vsite, constr,
1346 nrnb, wcycle);
1349 /* Take a trial step to this new point - new coords in s_b */
1350 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, &s_min->s.cg_p, s_b,
1351 constr, top, nrnb, wcycle, -1);
1353 neval++;
1354 /* Calculate energy for the trial step */
1355 evaluate_energy(fplog, cr,
1356 top_global, s_b, top,
1357 inputrec, nrnb, wcycle, gstat,
1358 vsite, constr, fcd, graph, mdatoms, fr,
1359 mu_tot, enerd, vir, pres, -1, FALSE);
1361 /* p does not change within a step, but since the domain decomposition
1362 * might change, we have to use cg_p of s_b here.
1364 const rvec *pb = as_rvec_array(s_b->s.cg_p.data());
1365 const rvec *sfb = as_rvec_array(s_b->f.data());
1366 gpb = 0;
1367 for (int i = 0; i < mdatoms->homenr; i++)
1369 for (m = 0; m < DIM; m++)
1371 gpb -= pb[i][m]*sfb[i][m]; /* f is negative gradient, thus the sign */
1374 /* Sum the gradient along the line across CPUs */
1375 if (PAR(cr))
1377 gmx_sumd(1, &gpb, cr);
1380 if (debug)
1382 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1383 s_a->epot, s_b->epot, s_c->epot, gpb);
1386 epot_repl = s_b->epot;
1388 /* Keep one of the intervals based on the value of the derivative at the new point */
1389 if (gpb > 0)
1391 /* Replace c endpoint with b */
1392 swap_em_state(&s_b, &s_c);
1393 c = b;
1394 gpc = gpb;
1396 else
1398 /* Replace a endpoint with b */
1399 swap_em_state(&s_b, &s_a);
1400 a = b;
1401 gpa = gpb;
1405 * Stop search as soon as we find a value smaller than the endpoints.
1406 * Never run more than 20 steps, no matter what.
1408 nminstep++;
1410 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1411 (nminstep < 20));
1413 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1414 nminstep >= 20)
1416 /* OK. We couldn't find a significantly lower energy.
1417 * If beta==0 this was steepest descent, and then we give up.
1418 * If not, set beta=0 and restart with steepest descent before quitting.
1420 if (beta == 0.0)
1422 /* Converged */
1423 converged = TRUE;
1424 break;
1426 else
1428 /* Reset memory before giving up */
1429 beta = 0.0;
1430 continue;
1434 /* Select min energy state of A & C, put the best in B.
1436 if (s_c->epot < s_a->epot)
1438 if (debug)
1440 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1441 s_c->epot, s_a->epot);
1443 swap_em_state(&s_b, &s_c);
1444 gpb = gpc;
1446 else
1448 if (debug)
1450 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1451 s_a->epot, s_c->epot);
1453 swap_em_state(&s_b, &s_a);
1454 gpb = gpa;
1458 else
1460 if (debug)
1462 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1463 s_c->epot);
1465 swap_em_state(&s_b, &s_c);
1466 gpb = gpc;
1469 /* new search direction */
1470 /* beta = 0 means forget all memory and restart with steepest descents. */
1471 if (nstcg && ((step % nstcg) == 0))
1473 beta = 0.0;
1475 else
1477 /* s_min->fnorm cannot be zero, because then we would have converged
1478 * and broken out.
1481 /* Polak-Ribiere update.
1482 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1484 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1486 /* Limit beta to prevent oscillations */
1487 if (fabs(beta) > 5.0)
1489 beta = 0.0;
1493 /* update positions */
1494 swap_em_state(&s_min, &s_b);
1495 gpa = gpb;
1497 /* Print it if necessary */
1498 if (MASTER(cr))
1500 if (bVerbose)
1502 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1503 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1504 step, s_min->epot, s_min->fnorm/sqrtNumAtoms,
1505 s_min->fmax, s_min->a_fmax+1);
1506 fflush(stderr);
1508 /* Store the new (lower) energies */
1509 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1510 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1511 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1513 do_log = do_per_step(step, inputrec->nstlog);
1514 do_ene = do_per_step(step, inputrec->nstenergy);
1516 /* Prepare IMD energy record, if bIMD is TRUE. */
1517 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1519 if (do_log)
1521 print_ebin_header(fplog, step, step);
1523 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1524 do_log ? fplog : nullptr, step, step, eprNORMAL,
1525 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1528 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1529 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, as_rvec_array(state_global->x.data()), inputrec, 0, wcycle) && MASTER(cr))
1531 IMD_send_positions(inputrec->imd);
1534 /* Stop when the maximum force lies below tolerance.
1535 * If we have reached machine precision, converged is already set to true.
1537 converged = converged || (s_min->fmax < inputrec->em_tol);
1539 } /* End of the loop */
1541 /* IMD cleanup, if bIMD is TRUE. */
1542 IMD_finalize(inputrec->bIMD, inputrec->imd);
1544 if (converged)
1546 step--; /* we never took that last step in this case */
1549 if (s_min->fmax > inputrec->em_tol)
1551 if (MASTER(cr))
1553 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1554 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1556 converged = FALSE;
1559 if (MASTER(cr))
1561 /* If we printed energy and/or logfile last step (which was the last step)
1562 * we don't have to do it again, but otherwise print the final values.
1564 if (!do_log)
1566 /* Write final value to log since we didn't do anything the last step */
1567 print_ebin_header(fplog, step, step);
1569 if (!do_ene || !do_log)
1571 /* Write final energy file entries */
1572 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1573 !do_log ? fplog : nullptr, step, step, eprNORMAL,
1574 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1578 /* Print some stuff... */
1579 if (MASTER(cr))
1581 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1584 /* IMPORTANT!
1585 * For accurate normal mode calculation it is imperative that we
1586 * store the last conformation into the full precision binary trajectory.
1588 * However, we should only do it if we did NOT already write this step
1589 * above (which we did if do_x or do_f was true).
1591 do_x = !do_per_step(step, inputrec->nstxout);
1592 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1594 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1595 top_global, inputrec, step,
1596 s_min, state_global, observablesHistory);
1599 if (MASTER(cr))
1601 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1602 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1603 s_min, sqrtNumAtoms);
1604 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1605 s_min, sqrtNumAtoms);
1607 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1610 finish_em(cr, outf, walltime_accounting, wcycle);
1612 /* To print the actual number of steps we needed somewhere */
1613 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1615 return 0;
1616 } /* That's all folks */
1619 /*! \brief Do L-BFGS conjugate gradients minimization
1620 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
1621 int nfile, const t_filenm fnm[],
1622 const gmx_output_env_t *oenv, gmx_bool bVerbose,
1623 int nstglobalcomm,
1624 gmx_vsite_t *vsite, gmx_constr_t constr,
1625 int stepout,
1626 gmx::IMDOutputProvider *outputProvider,
1627 t_inputrec *inputrec,
1628 gmx_mtop_t *top_global, t_fcdata *fcd,
1629 t_state *state_global,
1630 t_mdatoms *mdatoms,
1631 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1632 gmx_edsam_t ed,
1633 t_forcerec *fr,
1634 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
1635 gmx_membed_t gmx_unused *membed,
1636 real cpt_period, real max_hours,
1637 int imdport,
1638 unsigned long Flags,
1639 gmx_walltime_accounting_t walltime_accounting)
1641 double do_lbfgs(FILE *fplog, t_commrec *cr, const gmx::MDLogger gmx_unused &mdlog,
1642 int nfile, const t_filenm fnm[],
1643 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
1644 int gmx_unused nstglobalcomm,
1645 gmx_vsite_t *vsite, gmx_constr_t constr,
1646 int gmx_unused stepout,
1647 gmx::IMDOutputProvider *outputProvider,
1648 t_inputrec *inputrec,
1649 gmx_mtop_t *top_global, t_fcdata *fcd,
1650 t_state *state_global,
1651 ObservablesHistory *observablesHistory,
1652 t_mdatoms *mdatoms,
1653 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1654 gmx_edsam_t gmx_unused ed,
1655 t_forcerec *fr,
1656 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1657 gmx_membed_t gmx_unused *membed,
1658 real gmx_unused cpt_period, real gmx_unused max_hours,
1659 int imdport,
1660 unsigned long gmx_unused Flags,
1661 gmx_walltime_accounting_t walltime_accounting)
1663 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1664 em_state_t ems;
1665 gmx_localtop_t *top;
1666 gmx_enerdata_t *enerd;
1667 gmx_global_stat_t gstat;
1668 t_graph *graph;
1669 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1670 double stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1671 real *rho, *alpha, *p, *s, **dx, **dg;
1672 real a, b, c, maxdelta, delta;
1673 real diag, Epot0;
1674 real dgdx, dgdg, sq, yr, beta;
1675 t_mdebin *mdebin;
1676 gmx_bool converged;
1677 rvec mu_tot;
1678 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1679 tensor vir, pres;
1680 int start, end, number_steps;
1681 gmx_mdoutf_t outf;
1682 int i, k, m, n, gf, step;
1683 int mdof_flags;
1685 if (PAR(cr))
1687 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1690 if (nullptr != constr)
1692 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).");
1695 n = 3*state_global->natoms;
1696 nmaxcorr = inputrec->nbfgscorr;
1698 snew(frozen, n);
1700 snew(p, n);
1701 snew(rho, nmaxcorr);
1702 snew(alpha, nmaxcorr);
1704 snew(dx, nmaxcorr);
1705 for (i = 0; i < nmaxcorr; i++)
1707 snew(dx[i], n);
1710 snew(dg, nmaxcorr);
1711 for (i = 0; i < nmaxcorr; i++)
1713 snew(dg[i], n);
1716 step = 0;
1717 neval = 0;
1719 /* Init em */
1720 init_em(fplog, LBFGS, cr, outputProvider, inputrec,
1721 state_global, top_global, &ems, &top,
1722 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat,
1723 vsite, constr, nullptr,
1724 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
1726 start = 0;
1727 end = mdatoms->homenr;
1729 /* We need 4 working states */
1730 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1731 em_state_t *sa = &s0;
1732 em_state_t *sb = &s1;
1733 em_state_t *sc = &s2;
1734 em_state_t *last = &s3;
1735 /* Initialize by copying the state from ems (we could skip x and f here) */
1736 *sa = ems;
1737 *sb = ems;
1738 *sc = ems;
1740 /* Print to log file */
1741 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1743 do_log = do_ene = do_x = do_f = TRUE;
1745 /* Max number of steps */
1746 number_steps = inputrec->nsteps;
1748 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1749 gf = 0;
1750 for (i = start; i < end; i++)
1752 if (mdatoms->cFREEZE)
1754 gf = mdatoms->cFREEZE[i];
1756 for (m = 0; m < DIM; m++)
1758 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1761 if (MASTER(cr))
1763 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1765 if (fplog)
1767 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1770 if (vsite)
1772 construct_vsites(vsite, as_rvec_array(state_global->x.data()), 1, nullptr,
1773 top->idef.iparams, top->idef.il,
1774 fr->ePBC, fr->bMolPBC, cr, state_global->box);
1777 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1778 /* do_force always puts the charge groups in the box and shifts again
1779 * We do not unshift, so molecules are always whole
1781 neval++;
1782 evaluate_energy(fplog, cr,
1783 top_global, &ems, top,
1784 inputrec, nrnb, wcycle, gstat,
1785 vsite, constr, fcd, graph, mdatoms, fr,
1786 mu_tot, enerd, vir, pres, -1, TRUE);
1787 where();
1789 if (MASTER(cr))
1791 /* Copy stuff to the energy bin for easy printing etc. */
1792 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1793 mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
1794 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1796 print_ebin_header(fplog, step, step);
1797 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1798 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1800 where();
1802 /* Set the initial step.
1803 * since it will be multiplied by the non-normalized search direction
1804 * vector (force vector the first time), we scale it by the
1805 * norm of the force.
1808 if (MASTER(cr))
1810 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1811 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1812 fprintf(stderr, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1813 fprintf(stderr, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1814 fprintf(stderr, "\n");
1815 /* and copy to the log file too... */
1816 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1817 fprintf(fplog, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1818 fprintf(fplog, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1819 fprintf(fplog, "\n");
1822 // Point is an index to the memory of search directions, where 0 is the first one.
1823 point = 0;
1825 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1826 real *fInit = static_cast<real *>(as_rvec_array(ems.f.data())[0]);
1827 for (i = 0; i < n; i++)
1829 if (!frozen[i])
1831 dx[point][i] = fInit[i]; /* Initial search direction */
1833 else
1835 dx[point][i] = 0;
1839 // Stepsize will be modified during the search, and actually it is not critical
1840 // (the main efficiency in the algorithm comes from changing directions), but
1841 // we still need an initial value, so estimate it as the inverse of the norm
1842 // so we take small steps where the potential fluctuates a lot.
1843 stepsize = 1.0/ems.fnorm;
1845 /* Start the loop over BFGS steps.
1846 * Each successful step is counted, and we continue until
1847 * we either converge or reach the max number of steps.
1850 ncorr = 0;
1852 /* Set the gradient from the force */
1853 converged = FALSE;
1854 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1857 /* Write coordinates if necessary */
1858 do_x = do_per_step(step, inputrec->nstxout);
1859 do_f = do_per_step(step, inputrec->nstfout);
1861 mdof_flags = 0;
1862 if (do_x)
1864 mdof_flags |= MDOF_X;
1867 if (do_f)
1869 mdof_flags |= MDOF_F;
1872 if (inputrec->bIMD)
1874 mdof_flags |= MDOF_IMD;
1877 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1878 top_global, step, (real)step, &ems.s, state_global, observablesHistory, &ems.f);
1880 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1882 /* make s a pointer to current search direction - point=0 first time we get here */
1883 s = dx[point];
1885 real *xx = static_cast<real *>(as_rvec_array(ems.s.x.data())[0]);
1886 real *ff = static_cast<real *>(as_rvec_array(ems.f.data())[0]);
1888 // calculate line gradient in position A
1889 for (gpa = 0, i = 0; i < n; i++)
1891 gpa -= s[i]*ff[i];
1894 /* Calculate minimum allowed stepsize along the line, before the average (norm)
1895 * relative change in coordinate is smaller than precision
1897 for (minstep = 0, i = 0; i < n; i++)
1899 tmp = fabs(xx[i]);
1900 if (tmp < 1.0)
1902 tmp = 1.0;
1904 tmp = s[i]/tmp;
1905 minstep += tmp*tmp;
1907 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1909 if (stepsize < minstep)
1911 converged = TRUE;
1912 break;
1915 // Before taking any steps along the line, store the old position
1916 *last = ems;
1917 real *lastx = static_cast<real *>(as_rvec_array(last->s.x.data())[0]);
1918 real *lastf = static_cast<real *>(as_rvec_array(last->f.data())[0]);
1919 Epot0 = ems.epot;
1921 *sa = ems;
1923 /* Take a step downhill.
1924 * In theory, we should find the actual minimum of the function in this
1925 * direction, somewhere along the line.
1926 * That is quite possible, but it turns out to take 5-10 function evaluations
1927 * for each line. However, we dont really need to find the exact minimum -
1928 * it is much better to start a new BFGS step in a modified direction as soon
1929 * as we are close to it. This will save a lot of energy evaluations.
1931 * In practice, we just try to take a single step.
1932 * If it worked (i.e. lowered the energy), we increase the stepsize but
1933 * continue straight to the next BFGS step without trying to find any minimum,
1934 * i.e. we change the search direction too. If the line was smooth, it is
1935 * likely we are in a smooth region, and then it makes sense to take longer
1936 * steps in the modified search direction too.
1938 * If it didn't work (higher energy), there must be a minimum somewhere between
1939 * the old position and the new one. Then we need to start by finding a lower
1940 * value before we change search direction. Since the energy was apparently
1941 * quite rough, we need to decrease the step size.
1943 * Due to the finite numerical accuracy, it turns out that it is a good idea
1944 * to accept a SMALL increase in energy, if the derivative is still downhill.
1945 * This leads to lower final energies in the tests I've done. / Erik
1948 // State "A" is the first position along the line.
1949 // reference position along line is initially zero
1950 a = 0.0;
1952 // Check stepsize first. We do not allow displacements
1953 // larger than emstep.
1957 // Pick a new position C by adding stepsize to A.
1958 c = a + stepsize;
1960 // Calculate what the largest change in any individual coordinate
1961 // would be (translation along line * gradient along line)
1962 maxdelta = 0;
1963 for (i = 0; i < n; i++)
1965 delta = c*s[i];
1966 if (delta > maxdelta)
1968 maxdelta = delta;
1971 // If any displacement is larger than the stepsize limit, reduce the step
1972 if (maxdelta > inputrec->em_stepsize)
1974 stepsize *= 0.1;
1977 while (maxdelta > inputrec->em_stepsize);
1979 // Take a trial step and move the coordinate array xc[] to position C
1980 real *xc = static_cast<real *>(as_rvec_array(sc->s.x.data())[0]);
1981 for (i = 0; i < n; i++)
1983 xc[i] = lastx[i] + c*s[i];
1986 neval++;
1987 // Calculate energy for the trial step in position C
1988 evaluate_energy(fplog, cr,
1989 top_global, sc, top,
1990 inputrec, nrnb, wcycle, gstat,
1991 vsite, constr, fcd, graph, mdatoms, fr,
1992 mu_tot, enerd, vir, pres, step, FALSE);
1994 // Calc line gradient in position C
1995 real *fc = static_cast<real *>(as_rvec_array(sc->f.data())[0]);
1996 for (gpc = 0, i = 0; i < n; i++)
1998 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
2000 /* Sum the gradient along the line across CPUs */
2001 if (PAR(cr))
2003 gmx_sumd(1, &gpc, cr);
2006 // This is the max amount of increase in energy we tolerate.
2007 // By allowing VERY small changes (close to numerical precision) we
2008 // frequently find even better (lower) final energies.
2009 tmp = sqrt(GMX_REAL_EPS)*fabs(sa->epot);
2011 // Accept the step if the energy is lower in the new position C (compared to A),
2012 // or if it is not significantly higher and the line derivative is still negative.
2013 if (sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp)))
2015 // Great, we found a better energy. We no longer try to alter the
2016 // stepsize, but simply accept this new better position. The we select a new
2017 // search direction instead, which will be much more efficient than continuing
2018 // to take smaller steps along a line. Set fnorm based on the new C position,
2019 // which will be used to update the stepsize to 1/fnorm further down.
2020 foundlower = TRUE;
2022 else
2024 // If we got here, the energy is NOT lower in point C, i.e. it will be the same
2025 // or higher than in point A. In this case it is pointless to move to point C,
2026 // so we will have to do more iterations along the same line to find a smaller
2027 // value in the interval [A=0.0,C].
2028 // Here, A is still 0.0, but that will change when we do a search in the interval
2029 // [0.0,C] below. That search we will do by interpolation or bisection rather
2030 // than with the stepsize, so no need to modify it. For the next search direction
2031 // it will be reset to 1/fnorm anyway.
2032 foundlower = FALSE;
2035 if (!foundlower)
2037 // OK, if we didn't find a lower value we will have to locate one now - there must
2038 // be one in the interval [a,c].
2039 // The same thing is valid here, though: Don't spend dozens of iterations to find
2040 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2041 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2042 // I also have a safeguard for potentially really pathological functions so we never
2043 // take more than 20 steps before we give up.
2044 // If we already found a lower value we just skip this step and continue to the update.
2045 real fnorm = 0;
2046 nminstep = 0;
2049 // Select a new trial point B in the interval [A,C].
2050 // If the derivatives at points a & c have different sign we interpolate to zero,
2051 // otherwise just do a bisection since there might be multiple minima/maxima
2052 // inside the interval.
2053 if (gpa < 0 && gpc > 0)
2055 b = a + gpa*(a-c)/(gpc-gpa);
2057 else
2059 b = 0.5*(a+c);
2062 /* safeguard if interpolation close to machine accuracy causes errors:
2063 * never go outside the interval
2065 if (b <= a || b >= c)
2067 b = 0.5*(a+c);
2070 // Take a trial step to point B
2071 real *xb = static_cast<real *>(as_rvec_array(sb->s.x.data())[0]);
2072 for (i = 0; i < n; i++)
2074 xb[i] = lastx[i] + b*s[i];
2077 neval++;
2078 // Calculate energy for the trial step in point B
2079 evaluate_energy(fplog, cr,
2080 top_global, sb, top,
2081 inputrec, nrnb, wcycle, gstat,
2082 vsite, constr, fcd, graph, mdatoms, fr,
2083 mu_tot, enerd, vir, pres, step, FALSE);
2084 fnorm = sb->fnorm;
2086 // Calculate gradient in point B
2087 real *fb = static_cast<real *>(as_rvec_array(sb->f.data())[0]);
2088 for (gpb = 0, i = 0; i < n; i++)
2090 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2093 /* Sum the gradient along the line across CPUs */
2094 if (PAR(cr))
2096 gmx_sumd(1, &gpb, cr);
2099 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2100 // at the new point B, and rename the endpoints of this new interval A and C.
2101 if (gpb > 0)
2103 /* Replace c endpoint with b */
2104 c = b;
2105 /* swap states b and c */
2106 swap_em_state(&sb, &sc);
2108 else
2110 /* Replace a endpoint with b */
2111 a = b;
2112 /* swap states a and b */
2113 swap_em_state(&sa, &sb);
2117 * Stop search as soon as we find a value smaller than the endpoints,
2118 * or if the tolerance is below machine precision.
2119 * Never run more than 20 steps, no matter what.
2121 nminstep++;
2123 while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2125 if (fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2127 /* OK. We couldn't find a significantly lower energy.
2128 * If ncorr==0 this was steepest descent, and then we give up.
2129 * If not, reset memory to restart as steepest descent before quitting.
2131 if (ncorr == 0)
2133 /* Converged */
2134 converged = TRUE;
2135 break;
2137 else
2139 /* Reset memory */
2140 ncorr = 0;
2141 /* Search in gradient direction */
2142 for (i = 0; i < n; i++)
2144 dx[point][i] = ff[i];
2146 /* Reset stepsize */
2147 stepsize = 1.0/fnorm;
2148 continue;
2152 /* Select min energy state of A & C, put the best in xx/ff/Epot
2154 if (sc->epot < sa->epot)
2156 /* Use state C */
2157 ems = *sc;
2158 step_taken = c;
2160 else
2162 /* Use state A */
2163 ems = *sa;
2164 step_taken = a;
2168 else
2170 /* found lower */
2171 /* Use state C */
2172 ems = *sc;
2173 step_taken = c;
2176 /* Update the memory information, and calculate a new
2177 * approximation of the inverse hessian
2180 /* Have new data in Epot, xx, ff */
2181 if (ncorr < nmaxcorr)
2183 ncorr++;
2186 for (i = 0; i < n; i++)
2188 dg[point][i] = lastf[i]-ff[i];
2189 dx[point][i] *= step_taken;
2192 dgdg = 0;
2193 dgdx = 0;
2194 for (i = 0; i < n; i++)
2196 dgdg += dg[point][i]*dg[point][i];
2197 dgdx += dg[point][i]*dx[point][i];
2200 diag = dgdx/dgdg;
2202 rho[point] = 1.0/dgdx;
2203 point++;
2205 if (point >= nmaxcorr)
2207 point = 0;
2210 /* Update */
2211 for (i = 0; i < n; i++)
2213 p[i] = ff[i];
2216 cp = point;
2218 /* Recursive update. First go back over the memory points */
2219 for (k = 0; k < ncorr; k++)
2221 cp--;
2222 if (cp < 0)
2224 cp = ncorr-1;
2227 sq = 0;
2228 for (i = 0; i < n; i++)
2230 sq += dx[cp][i]*p[i];
2233 alpha[cp] = rho[cp]*sq;
2235 for (i = 0; i < n; i++)
2237 p[i] -= alpha[cp]*dg[cp][i];
2241 for (i = 0; i < n; i++)
2243 p[i] *= diag;
2246 /* And then go forward again */
2247 for (k = 0; k < ncorr; k++)
2249 yr = 0;
2250 for (i = 0; i < n; i++)
2252 yr += p[i]*dg[cp][i];
2255 beta = rho[cp]*yr;
2256 beta = alpha[cp]-beta;
2258 for (i = 0; i < n; i++)
2260 p[i] += beta*dx[cp][i];
2263 cp++;
2264 if (cp >= ncorr)
2266 cp = 0;
2270 for (i = 0; i < n; i++)
2272 if (!frozen[i])
2274 dx[point][i] = p[i];
2276 else
2278 dx[point][i] = 0;
2282 /* Print it if necessary */
2283 if (MASTER(cr))
2285 if (bVerbose)
2287 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2288 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2289 step, ems.epot, ems.fnorm/sqrtNumAtoms, ems.fmax, ems.a_fmax + 1);
2290 fflush(stderr);
2292 /* Store the new (lower) energies */
2293 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2294 mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
2295 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2296 do_log = do_per_step(step, inputrec->nstlog);
2297 do_ene = do_per_step(step, inputrec->nstenergy);
2298 if (do_log)
2300 print_ebin_header(fplog, step, step);
2302 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2303 do_log ? fplog : nullptr, step, step, eprNORMAL,
2304 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2307 /* Send x and E to IMD client, if bIMD is TRUE. */
2308 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, as_rvec_array(state_global->x.data()), inputrec, 0, wcycle) && MASTER(cr))
2310 IMD_send_positions(inputrec->imd);
2313 // Reset stepsize in we are doing more iterations
2314 stepsize = 1.0/ems.fnorm;
2316 /* Stop when the maximum force lies below tolerance.
2317 * If we have reached machine precision, converged is already set to true.
2319 converged = converged || (ems.fmax < inputrec->em_tol);
2321 } /* End of the loop */
2323 /* IMD cleanup, if bIMD is TRUE. */
2324 IMD_finalize(inputrec->bIMD, inputrec->imd);
2326 if (converged)
2328 step--; /* we never took that last step in this case */
2331 if (ems.fmax > inputrec->em_tol)
2333 if (MASTER(cr))
2335 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2336 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2338 converged = FALSE;
2341 /* If we printed energy and/or logfile last step (which was the last step)
2342 * we don't have to do it again, but otherwise print the final values.
2344 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2346 print_ebin_header(fplog, step, step);
2348 if (!do_ene || !do_log) /* Write final energy file entries */
2350 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2351 !do_log ? fplog : nullptr, step, step, eprNORMAL,
2352 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2355 /* Print some stuff... */
2356 if (MASTER(cr))
2358 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2361 /* IMPORTANT!
2362 * For accurate normal mode calculation it is imperative that we
2363 * store the last conformation into the full precision binary trajectory.
2365 * However, we should only do it if we did NOT already write this step
2366 * above (which we did if do_x or do_f was true).
2368 do_x = !do_per_step(step, inputrec->nstxout);
2369 do_f = !do_per_step(step, inputrec->nstfout);
2370 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2371 top_global, inputrec, step,
2372 &ems, state_global, observablesHistory);
2374 if (MASTER(cr))
2376 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2377 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2378 number_steps, &ems, sqrtNumAtoms);
2379 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2380 number_steps, &ems, sqrtNumAtoms);
2382 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2385 finish_em(cr, outf, walltime_accounting, wcycle);
2387 /* To print the actual number of steps we needed somewhere */
2388 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2390 return 0;
2391 } /* That's all folks */
2393 /*! \brief Do steepest descents minimization
2394 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
2395 int nfile, const t_filenm fnm[],
2396 const gmx_output_env_t *oenv, gmx_bool bVerbose,
2397 int nstglobalcomm,
2398 gmx_vsite_t *vsite, gmx_constr_t constr,
2399 int stepout,
2400 gmx::IMDOutputProvider *outputProvider,
2401 t_inputrec *inputrec,
2402 gmx_mtop_t *top_global, t_fcdata *fcd,
2403 t_state *state_global,
2404 t_mdatoms *mdatoms,
2405 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2406 gmx_edsam_t ed,
2407 t_forcerec *fr,
2408 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
2409 real cpt_period, real max_hours,
2410 int imdport,
2411 unsigned long Flags,
2412 gmx_walltime_accounting_t walltime_accounting)
2414 double do_steep(FILE *fplog, t_commrec *cr, const gmx::MDLogger gmx_unused &mdlog,
2415 int nfile, const t_filenm fnm[],
2416 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
2417 int gmx_unused nstglobalcomm,
2418 gmx_vsite_t *vsite, gmx_constr_t constr,
2419 int gmx_unused stepout,
2420 gmx::IMDOutputProvider *outputProvider,
2421 t_inputrec *inputrec,
2422 gmx_mtop_t *top_global, t_fcdata *fcd,
2423 t_state *state_global,
2424 ObservablesHistory *observablesHistory,
2425 t_mdatoms *mdatoms,
2426 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2427 gmx_edsam_t gmx_unused ed,
2428 t_forcerec *fr,
2429 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2430 gmx_membed_t gmx_unused *membed,
2431 real gmx_unused cpt_period, real gmx_unused max_hours,
2432 int imdport,
2433 unsigned long gmx_unused Flags,
2434 gmx_walltime_accounting_t walltime_accounting)
2436 const char *SD = "Steepest Descents";
2437 gmx_localtop_t *top;
2438 gmx_enerdata_t *enerd;
2439 gmx_global_stat_t gstat;
2440 t_graph *graph;
2441 real stepsize;
2442 real ustep;
2443 gmx_mdoutf_t outf;
2444 t_mdebin *mdebin;
2445 gmx_bool bDone, bAbort, do_x, do_f;
2446 tensor vir, pres;
2447 rvec mu_tot;
2448 int nsteps;
2449 int count = 0;
2450 int steps_accepted = 0;
2452 /* Create 2 states on the stack and extract pointers that we will swap */
2453 em_state_t s0 {}, s1 {};
2454 em_state_t *s_min = &s0;
2455 em_state_t *s_try = &s1;
2457 /* Init em and store the local state in s_try */
2458 init_em(fplog, SD, cr, outputProvider, inputrec,
2459 state_global, top_global, s_try, &top,
2460 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat,
2461 vsite, constr, nullptr,
2462 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
2464 /* Print to log file */
2465 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2467 /* Set variables for stepsize (in nm). This is the largest
2468 * step that we are going to make in any direction.
2470 ustep = inputrec->em_stepsize;
2471 stepsize = 0;
2473 /* Max number of steps */
2474 nsteps = inputrec->nsteps;
2476 if (MASTER(cr))
2478 /* Print to the screen */
2479 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2481 if (fplog)
2483 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2486 /**** HERE STARTS THE LOOP ****
2487 * count is the counter for the number of steps
2488 * bDone will be TRUE when the minimization has converged
2489 * bAbort will be TRUE when nsteps steps have been performed or when
2490 * the stepsize becomes smaller than is reasonable for machine precision
2492 count = 0;
2493 bDone = FALSE;
2494 bAbort = FALSE;
2495 while (!bDone && !bAbort)
2497 bAbort = (nsteps >= 0) && (count == nsteps);
2499 /* set new coordinates, except for first step */
2500 bool validStep = true;
2501 if (count > 0)
2503 validStep =
2504 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2505 s_min, stepsize, &s_min->f, s_try,
2506 constr, top, nrnb, wcycle, count);
2509 if (validStep)
2511 evaluate_energy(fplog, cr,
2512 top_global, s_try, top,
2513 inputrec, nrnb, wcycle, gstat,
2514 vsite, constr, fcd, graph, mdatoms, fr,
2515 mu_tot, enerd, vir, pres, count, count == 0);
2517 else
2519 // Signal constraint error during stepping with energy=inf
2520 s_try->epot = std::numeric_limits<real>::infinity();
2523 if (MASTER(cr))
2525 print_ebin_header(fplog, count, count);
2528 if (count == 0)
2530 s_min->epot = s_try->epot;
2533 /* Print it if necessary */
2534 if (MASTER(cr))
2536 if (bVerbose)
2538 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2539 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2540 ( (count == 0) || (s_try->epot < s_min->epot) ) ? '\n' : '\r');
2541 fflush(stderr);
2544 if ( (count == 0) || (s_try->epot < s_min->epot) )
2546 /* Store the new (lower) energies */
2547 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2548 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2549 s_try->s.box, nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2551 /* Prepare IMD energy record, if bIMD is TRUE. */
2552 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2554 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2555 do_per_step(steps_accepted, inputrec->nstdisreout),
2556 do_per_step(steps_accepted, inputrec->nstorireout),
2557 fplog, count, count, eprNORMAL,
2558 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2559 fflush(fplog);
2563 /* Now if the new energy is smaller than the previous...
2564 * or if this is the first step!
2565 * or if we did random steps!
2568 if ( (count == 0) || (s_try->epot < s_min->epot) )
2570 steps_accepted++;
2572 /* Test whether the convergence criterion is met... */
2573 bDone = (s_try->fmax < inputrec->em_tol);
2575 /* Copy the arrays for force, positions and energy */
2576 /* The 'Min' array always holds the coords and forces of the minimal
2577 sampled energy */
2578 swap_em_state(&s_min, &s_try);
2579 if (count > 0)
2581 ustep *= 1.2;
2584 /* Write to trn, if necessary */
2585 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2586 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2587 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
2588 top_global, inputrec, count,
2589 s_min, state_global, observablesHistory);
2591 else
2593 /* If energy is not smaller make the step smaller... */
2594 ustep *= 0.5;
2596 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2598 /* Reload the old state */
2599 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2600 s_min, top, mdatoms, fr, vsite, constr,
2601 nrnb, wcycle);
2605 /* Determine new step */
2606 stepsize = ustep/s_min->fmax;
2608 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2609 #if GMX_DOUBLE
2610 if (count == nsteps || ustep < 1e-12)
2611 #else
2612 if (count == nsteps || ustep < 1e-6)
2613 #endif
2615 if (MASTER(cr))
2617 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != nullptr);
2618 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != nullptr);
2620 bAbort = TRUE;
2623 /* Send IMD energies and positions, if bIMD is TRUE. */
2624 if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box, as_rvec_array(state_global->x.data()), inputrec, 0, wcycle) && MASTER(cr))
2626 IMD_send_positions(inputrec->imd);
2629 count++;
2630 } /* End of the loop */
2632 /* IMD cleanup, if bIMD is TRUE. */
2633 IMD_finalize(inputrec->bIMD, inputrec->imd);
2635 /* Print some data... */
2636 if (MASTER(cr))
2638 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2640 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2641 top_global, inputrec, count,
2642 s_min, state_global, observablesHistory);
2644 if (MASTER(cr))
2646 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2648 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2649 s_min, sqrtNumAtoms);
2650 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2651 s_min, sqrtNumAtoms);
2654 finish_em(cr, outf, walltime_accounting, wcycle);
2656 /* To print the actual number of steps we needed somewhere */
2657 inputrec->nsteps = count;
2659 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2661 return 0;
2662 } /* That's all folks */
2664 /*! \brief Do normal modes analysis
2665 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
2666 int nfile, const t_filenm fnm[],
2667 const gmx_output_env_t *oenv, gmx_bool bVerbose,
2668 int nstglobalcomm,
2669 gmx_vsite_t *vsite, gmx_constr_t constr,
2670 int stepout,
2671 gmx::IMDOutputProvider *outputProvider,
2672 t_inputrec *inputrec,
2673 gmx_mtop_t *top_global, t_fcdata *fcd,
2674 t_state *state_global,
2675 t_mdatoms *mdatoms,
2676 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2677 gmx_edsam_t ed,
2678 t_forcerec *fr,
2679 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
2680 real cpt_period, real max_hours,
2681 int imdport,
2682 unsigned long Flags,
2683 gmx_walltime_accounting_t walltime_accounting)
2685 double do_nm(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
2686 int nfile, const t_filenm fnm[],
2687 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
2688 int gmx_unused nstglobalcomm,
2689 gmx_vsite_t *vsite, gmx_constr_t constr,
2690 int gmx_unused stepout,
2691 gmx::IMDOutputProvider *outputProvider,
2692 t_inputrec *inputrec,
2693 gmx_mtop_t *top_global, t_fcdata *fcd,
2694 t_state *state_global,
2695 ObservablesHistory gmx_unused *observablesHistory,
2696 t_mdatoms *mdatoms,
2697 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2698 gmx_edsam_t gmx_unused ed,
2699 t_forcerec *fr,
2700 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2701 gmx_membed_t gmx_unused *membed,
2702 real gmx_unused cpt_period, real gmx_unused max_hours,
2703 int imdport,
2704 unsigned long gmx_unused Flags,
2705 gmx_walltime_accounting_t walltime_accounting)
2707 const char *NM = "Normal Mode Analysis";
2708 gmx_mdoutf_t outf;
2709 int nnodes, node;
2710 gmx_localtop_t *top;
2711 gmx_enerdata_t *enerd;
2712 gmx_global_stat_t gstat;
2713 t_graph *graph;
2714 tensor vir, pres;
2715 rvec mu_tot;
2716 rvec *fneg, *dfdx;
2717 gmx_bool bSparse; /* use sparse matrix storage format */
2718 size_t sz;
2719 gmx_sparsematrix_t * sparse_matrix = nullptr;
2720 real * full_matrix = nullptr;
2722 /* added with respect to mdrun */
2723 int row, col;
2724 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2725 real x_min;
2726 bool bIsMaster = MASTER(cr);
2728 if (constr != nullptr)
2730 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2733 gmx_shellfc_t *shellfc;
2735 em_state_t state_work {};
2737 /* Init em and store the local state in state_minimum */
2738 init_em(fplog, NM, cr, outputProvider, inputrec,
2739 state_global, top_global, &state_work, &top,
2740 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat,
2741 vsite, constr, &shellfc,
2742 nfile, fnm, &outf, nullptr, imdport, Flags, wcycle);
2744 std::vector<size_t> atom_index = get_atom_index(top_global);
2745 snew(fneg, atom_index.size());
2746 snew(dfdx, atom_index.size());
2748 #if !GMX_DOUBLE
2749 if (bIsMaster)
2751 fprintf(stderr,
2752 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2753 " which MIGHT not be accurate enough for normal mode analysis.\n"
2754 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
2755 " are fairly modest even if you recompile in double precision.\n\n");
2757 #endif
2759 /* Check if we can/should use sparse storage format.
2761 * Sparse format is only useful when the Hessian itself is sparse, which it
2762 * will be when we use a cutoff.
2763 * For small systems (n<1000) it is easier to always use full matrix format, though.
2765 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2767 GMX_LOG(mdlog.warning).appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2768 bSparse = FALSE;
2770 else if (atom_index.size() < 1000)
2772 GMX_LOG(mdlog.warning).appendTextFormatted("Small system size (N=%d), using full Hessian format.",
2773 atom_index.size());
2774 bSparse = FALSE;
2776 else
2778 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2779 bSparse = TRUE;
2782 /* Number of dimensions, based on real atoms, that is not vsites or shell */
2783 sz = DIM*atom_index.size();
2785 fprintf(stderr, "Allocating Hessian memory...\n\n");
2787 if (bSparse)
2789 sparse_matrix = gmx_sparsematrix_init(sz);
2790 sparse_matrix->compressed_symmetric = TRUE;
2792 else
2794 snew(full_matrix, sz*sz);
2797 init_nrnb(nrnb);
2799 where();
2801 /* Write start time and temperature */
2802 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2804 /* fudge nr of steps to nr of atoms */
2805 inputrec->nsteps = atom_index.size()*2;
2807 if (bIsMaster)
2809 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2810 *(top_global->name), (int)inputrec->nsteps);
2813 nnodes = cr->nnodes;
2815 /* Make evaluate_energy do a single node force calculation */
2816 cr->nnodes = 1;
2817 evaluate_energy(fplog, cr,
2818 top_global, &state_work, top,
2819 inputrec, nrnb, wcycle, gstat,
2820 vsite, constr, fcd, graph, mdatoms, fr,
2821 mu_tot, enerd, vir, pres, -1, TRUE);
2822 cr->nnodes = nnodes;
2824 /* if forces are not small, warn user */
2825 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
2827 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
2828 if (state_work.fmax > 1.0e-3)
2830 GMX_LOG(mdlog.warning).appendText(
2831 "The force is probably not small enough to "
2832 "ensure that you are at a minimum.\n"
2833 "Be aware that negative eigenvalues may occur\n"
2834 "when the resulting matrix is diagonalized.");
2837 /***********************************************************
2839 * Loop over all pairs in matrix
2841 * do_force called twice. Once with positive and
2842 * once with negative displacement
2844 ************************************************************/
2846 /* Steps are divided one by one over the nodes */
2847 bool bNS = true;
2848 for (unsigned int aid = cr->nodeid; aid < atom_index.size(); aid += nnodes)
2850 size_t atom = atom_index[aid];
2851 for (size_t d = 0; d < DIM; d++)
2853 gmx_bool bBornRadii = FALSE;
2854 gmx_int64_t step = 0;
2855 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2856 double t = 0;
2858 x_min = state_work.s.x[atom][d];
2860 for (unsigned int dx = 0; (dx < 2); dx++)
2862 if (dx == 0)
2864 state_work.s.x[atom][d] = x_min - der_range;
2866 else
2868 state_work.s.x[atom][d] = x_min + der_range;
2871 /* Make evaluate_energy do a single node force calculation */
2872 cr->nnodes = 1;
2873 if (shellfc)
2875 /* Now is the time to relax the shells */
2876 (void) relax_shell_flexcon(fplog, cr, bVerbose, step,
2877 inputrec, bNS, force_flags,
2878 top,
2879 constr, enerd, fcd,
2880 &state_work.s, &state_work.f, vir, mdatoms,
2881 nrnb, wcycle, graph, &top_global->groups,
2882 shellfc, fr, bBornRadii, t, mu_tot,
2883 vsite);
2884 bNS = false;
2885 step++;
2887 else
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, atom*2+dx, FALSE);
2896 cr->nnodes = nnodes;
2898 if (dx == 0)
2900 for (size_t i = 0; i < atom_index.size(); i++)
2902 copy_rvec(state_work.f[atom_index[i]], fneg[i]);
2907 /* x is restored to original */
2908 state_work.s.x[atom][d] = x_min;
2910 for (size_t j = 0; j < atom_index.size(); j++)
2912 for (size_t k = 0; (k < DIM); k++)
2914 dfdx[j][k] =
2915 -(state_work.f[atom_index[j]][k] - fneg[j][k])/(2*der_range);
2919 if (!bIsMaster)
2921 #if GMX_MPI
2922 #define mpi_type GMX_MPI_REAL
2923 MPI_Send(dfdx[0], atom_index.size()*DIM, mpi_type, MASTER(cr),
2924 cr->nodeid, cr->mpi_comm_mygroup);
2925 #endif
2927 else
2929 for (node = 0; (node < nnodes && atom+node < atom_index.size()); node++)
2931 if (node > 0)
2933 #if GMX_MPI
2934 MPI_Status stat;
2935 MPI_Recv(dfdx[0], atom_index.size()*DIM, mpi_type, node, node,
2936 cr->mpi_comm_mygroup, &stat);
2937 #undef mpi_type
2938 #endif
2941 row = (atom + node)*DIM + d;
2943 for (size_t j = 0; j < atom_index.size(); j++)
2945 for (size_t k = 0; k < DIM; k++)
2947 col = j*DIM + k;
2949 if (bSparse)
2951 if (col >= row && dfdx[j][k] != 0.0)
2953 gmx_sparsematrix_increment_value(sparse_matrix,
2954 row, col, dfdx[j][k]);
2957 else
2959 full_matrix[row*sz+col] = dfdx[j][k];
2966 if (bVerbose && fplog)
2968 fflush(fplog);
2971 /* write progress */
2972 if (bIsMaster && bVerbose)
2974 fprintf(stderr, "\rFinished step %d out of %d",
2975 static_cast<int>(std::min(atom+nnodes, atom_index.size())),
2976 static_cast<int>(atom_index.size()));
2977 fflush(stderr);
2981 if (bIsMaster)
2983 fprintf(stderr, "\n\nWriting Hessian...\n");
2984 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2987 finish_em(cr, outf, walltime_accounting, wcycle);
2989 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size()*2);
2991 return 0;
2994 } // namespace gmx