Improved state_change_natoms and used it more
[gromacs.git] / src / gromacs / mdlib / minimize.cpp
blobdd3af64376292b4b6178ac22f3202ee4f9889d48
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 void init_em(FILE *fplog, const char *title,
325 t_commrec *cr, t_inputrec *ir,
326 t_state *state_global, gmx_mtop_t *top_global,
327 em_state_t *ems, gmx_localtop_t **top,
328 t_nrnb *nrnb, rvec mu_tot,
329 t_forcerec *fr, gmx_enerdata_t **enerd,
330 t_graph **graph, t_mdatoms *mdatoms, gmx_global_stat_t *gstat,
331 gmx_vsite_t *vsite, gmx_constr_t constr, gmx_shellfc_t **shellfc,
332 int nfile, const t_filenm fnm[],
333 gmx_mdoutf_t *outf, t_mdebin **mdebin,
334 int imdport, unsigned long gmx_unused Flags,
335 gmx_wallcycle_t wcycle)
337 real dvdl_constr;
339 if (fplog)
341 fprintf(fplog, "Initiating %s\n", title);
344 state_global->ngtc = 0;
346 /* Initialize lambda variables */
347 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, nullptr);
349 init_nrnb(nrnb);
351 /* Interactive molecular dynamics */
352 init_IMD(ir, cr, top_global, fplog, 1, as_rvec_array(state_global->x.data()),
353 nfile, fnm, nullptr, imdport, Flags);
355 if (ir->eI == eiNM)
357 GMX_ASSERT(shellfc != NULL, "With NM we always support shells");
359 *shellfc = init_shell_flexcon(stdout,
360 top_global,
361 n_flexible_constraints(constr),
362 ir->nstcalcenergy,
363 DOMAINDECOMP(cr));
365 else
367 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI), "This else currently only handles energy minimizers, consider if your algorithm needs shell/flexible-constraint support");
369 /* With energy minimization, shells and flexible constraints are
370 * automatically minimized when treated like normal DOFS.
372 if (shellfc != nullptr)
374 *shellfc = nullptr;
378 if (DOMAINDECOMP(cr))
380 *top = dd_init_local_top(top_global);
382 dd_init_local_state(cr->dd, state_global, &ems->s);
384 /* Distribute the charge groups over the nodes from the master node */
385 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
386 state_global, top_global, ir,
387 &ems->s, &ems->f, mdatoms, *top,
388 fr, vsite, constr,
389 nrnb, nullptr, FALSE);
390 dd_store_state(cr->dd, &ems->s);
392 *graph = nullptr;
394 else
396 state_change_natoms(state_global, state_global->natoms);
397 /* Just copy the state */
398 ems->s = *state_global;
399 state_change_natoms(&ems->s, ems->s.natoms);
400 /* We need to allocate one element extra, since we might use
401 * (unaligned) 4-wide SIMD loads to access rvec entries.
403 ems->f.resize(ems->s.natoms + 1);
405 snew(*top, 1);
406 mdAlgorithmsSetupAtomData(cr, ir, top_global, *top, fr,
407 graph, mdatoms,
408 vsite, shellfc ? *shellfc : nullptr);
410 if (vsite)
412 set_vsite_top(vsite, *top, mdatoms, cr);
416 update_mdatoms(mdatoms, state_global->lambda[efptMASS]);
418 if (constr)
420 if (ir->eConstrAlg == econtSHAKE &&
421 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
423 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
424 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
427 if (!DOMAINDECOMP(cr))
429 set_constraints(constr, *top, ir, mdatoms, cr);
432 if (!ir->bContinuation)
434 /* Constrain the starting coordinates */
435 dvdl_constr = 0;
436 constrain(PAR(cr) ? nullptr : fplog, TRUE, TRUE, constr, &(*top)->idef,
437 ir, cr, -1, 0, 1.0, mdatoms,
438 as_rvec_array(ems->s.x.data()),
439 as_rvec_array(ems->s.x.data()),
440 nullptr,
441 fr->bMolPBC, ems->s.box,
442 ems->s.lambda[efptFEP], &dvdl_constr,
443 nullptr, nullptr, nrnb, econqCoord);
447 if (PAR(cr))
449 *gstat = global_stat_init(ir);
451 else
453 *gstat = nullptr;
456 *outf = init_mdoutf(fplog, nfile, fnm, 0, cr, ir, top_global, nullptr, wcycle);
458 snew(*enerd, 1);
459 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
460 *enerd);
462 if (mdebin != nullptr)
464 /* Init bin for energy stuff */
465 *mdebin = init_mdebin(mdoutf_get_fp_ene(*outf), top_global, ir, nullptr);
468 clear_rvec(mu_tot);
469 calc_shifts(ems->s.box, fr->shift_vec);
472 //! Finalize the minimization
473 static void finish_em(t_commrec *cr, gmx_mdoutf_t outf, const t_inputrec *ir,
474 gmx_walltime_accounting_t walltime_accounting,
475 gmx_wallcycle_t wcycle)
477 if (!(cr->duty & DUTY_PME))
479 /* Tell the PME only node to finish */
480 gmx_pme_send_finish(cr);
483 done_mdoutf(outf, ir);
485 em_time_end(walltime_accounting, wcycle);
488 //! Swap two different EM states during minimization
489 static void swap_em_state(em_state_t **ems1, em_state_t **ems2)
491 em_state_t *tmp;
493 tmp = *ems1;
494 *ems1 = *ems2;
495 *ems2 = tmp;
498 //! Save the EM trajectory
499 static void write_em_traj(FILE *fplog, t_commrec *cr,
500 gmx_mdoutf_t outf,
501 gmx_bool bX, gmx_bool bF, const char *confout,
502 gmx_mtop_t *top_global,
503 t_inputrec *ir, gmx_int64_t step,
504 em_state_t *state,
505 t_state *state_global,
506 energyhistory_t *energyHistory)
508 int mdof_flags = 0;
510 if (bX)
512 mdof_flags |= MDOF_X;
514 if (bF)
516 mdof_flags |= MDOF_F;
519 /* If we want IMD output, set appropriate MDOF flag */
520 if (ir->bIMD)
522 mdof_flags |= MDOF_IMD;
525 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
526 top_global, step, (double)step,
527 &state->s, state_global, energyHistory,
528 &state->f);
530 if (confout != nullptr && MASTER(cr))
532 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
534 /* Make molecules whole only for confout writing */
535 do_pbc_mtop(fplog, ir->ePBC, state_global->box, top_global,
536 as_rvec_array(state_global->x.data()));
539 write_sto_conf_mtop(confout,
540 *top_global->name, top_global,
541 as_rvec_array(state_global->x.data()), nullptr, ir->ePBC, state_global->box);
545 //! \brief Do one minimization step
547 // \returns true when the step succeeded, false when a constraint error occurred
548 static bool do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
549 gmx_bool bMolPBC,
550 em_state_t *ems1, real a, const PaddedRVecVector *force,
551 em_state_t *ems2,
552 gmx_constr_t constr, gmx_localtop_t *top,
553 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
554 gmx_int64_t count)
557 t_state *s1, *s2;
558 int start, end;
559 real dvdl_constr;
560 int nthreads gmx_unused;
562 bool validStep = true;
564 s1 = &ems1->s;
565 s2 = &ems2->s;
567 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
569 gmx_incons("state mismatch in do_em_step");
572 s2->flags = s1->flags;
574 if (s2->natoms != s1->natoms)
576 state_change_natoms(s2, s1->natoms);
577 /* We need to allocate one element extra, since we might use
578 * (unaligned) 4-wide SIMD loads to access rvec entries.
580 ems2->f.resize(s2->natoms + 1);
582 if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
584 s2->cg_gl.resize(s1->cg_gl.size());
587 copy_mat(s1->box, s2->box);
588 /* Copy free energy state */
589 s2->lambda = s1->lambda;
590 copy_mat(s1->box, s2->box);
592 start = 0;
593 end = md->homenr;
595 // cppcheck-suppress unreadVariable
596 nthreads = gmx_omp_nthreads_get(emntUpdate);
597 #pragma omp parallel num_threads(nthreads)
599 const rvec *x1 = as_rvec_array(s1->x.data());
600 rvec *x2 = as_rvec_array(s2->x.data());
601 const rvec *f = as_rvec_array(force->data());
603 int gf = 0;
604 #pragma omp for schedule(static) nowait
605 for (int i = start; i < end; i++)
609 if (md->cFREEZE)
611 gf = md->cFREEZE[i];
613 for (int m = 0; m < DIM; m++)
615 if (ir->opts.nFreeze[gf][m])
617 x2[i][m] = x1[i][m];
619 else
621 x2[i][m] = x1[i][m] + a*f[i][m];
625 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
628 if (s2->flags & (1<<estCGP))
630 /* Copy the CG p vector */
631 const rvec *p1 = as_rvec_array(s1->cg_p.data());
632 rvec *p2 = as_rvec_array(s2->cg_p.data());
633 #pragma omp for schedule(static) nowait
634 for (int i = start; i < end; i++)
636 // Trivial OpenMP block that does not throw
637 copy_rvec(p1[i], p2[i]);
641 if (DOMAINDECOMP(cr))
643 s2->ddp_count = s1->ddp_count;
645 /* OpenMP does not supported unsigned loop variables */
646 #pragma omp for schedule(static) nowait
647 for (int i = 0; i < static_cast<int>(s2->cg_gl.size()); i++)
649 s2->cg_gl[i] = s1->cg_gl[i];
651 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
655 if (constr)
657 wallcycle_start(wcycle, ewcCONSTR);
658 dvdl_constr = 0;
659 validStep =
660 constrain(nullptr, TRUE, TRUE, constr, &top->idef,
661 ir, cr, count, 0, 1.0, md,
662 as_rvec_array(s1->x.data()), as_rvec_array(s2->x.data()),
663 nullptr, bMolPBC, s2->box,
664 s2->lambda[efptBONDED], &dvdl_constr,
665 nullptr, nullptr, nrnb, econqCoord);
666 wallcycle_stop(wcycle, ewcCONSTR);
668 // We should move this check to the different minimizers
669 if (!validStep && ir->eI != eiSteep)
671 gmx_fatal(FARGS, "The coordinates could not be constrained. Minimizer '%s' can not handle constraint failures, use minimizer '%s' before using '%s'.",
672 EI(ir->eI), EI(eiSteep), EI(ir->eI));
676 return validStep;
679 //! Prepare EM for using domain decomposition parallellization
680 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
681 gmx_mtop_t *top_global, t_inputrec *ir,
682 em_state_t *ems, gmx_localtop_t *top,
683 t_mdatoms *mdatoms, t_forcerec *fr,
684 gmx_vsite_t *vsite, gmx_constr_t constr,
685 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
687 /* Repartition the domain decomposition */
688 dd_partition_system(fplog, step, cr, FALSE, 1,
689 nullptr, top_global, ir,
690 &ems->s, &ems->f,
691 mdatoms, top, fr, vsite, constr,
692 nrnb, wcycle, FALSE);
693 dd_store_state(cr->dd, &ems->s);
696 //! De one energy evaluation
697 static void evaluate_energy(FILE *fplog, t_commrec *cr,
698 gmx_mtop_t *top_global,
699 em_state_t *ems, gmx_localtop_t *top,
700 t_inputrec *inputrec,
701 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
702 gmx_global_stat_t gstat,
703 gmx_vsite_t *vsite, gmx_constr_t constr,
704 t_fcdata *fcd,
705 t_graph *graph, t_mdatoms *mdatoms,
706 t_forcerec *fr, rvec mu_tot,
707 gmx_enerdata_t *enerd, tensor vir, tensor pres,
708 gmx_int64_t count, gmx_bool bFirst)
710 real t;
711 gmx_bool bNS;
712 tensor force_vir, shake_vir, ekin;
713 real dvdl_constr, prescorr, enercorr, dvdlcorr;
714 real terminate = 0;
716 /* Set the time to the initial time, the time does not change during EM */
717 t = inputrec->init_t;
719 if (bFirst ||
720 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
722 /* This is the first state or an old state used before the last ns */
723 bNS = TRUE;
725 else
727 bNS = FALSE;
728 if (inputrec->nstlist > 0)
730 bNS = TRUE;
734 if (vsite)
736 construct_vsites(vsite, as_rvec_array(ems->s.x.data()), 1, nullptr,
737 top->idef.iparams, top->idef.il,
738 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
741 if (DOMAINDECOMP(cr) && bNS)
743 /* Repartition the domain decomposition */
744 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
745 ems, top, mdatoms, fr, vsite, constr,
746 nrnb, wcycle);
749 /* Calc force & energy on new trial position */
750 /* do_force always puts the charge groups in the box and shifts again
751 * We do not unshift, so molecules are always whole in congrad.c
753 do_force(fplog, cr, inputrec,
754 count, nrnb, wcycle, top, &top_global->groups,
755 ems->s.box, &ems->s.x, &ems->s.hist,
756 &ems->f, force_vir, mdatoms, enerd, fcd,
757 ems->s.lambda, graph, fr, vsite, mu_tot, t, nullptr, TRUE,
758 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
759 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
760 (bNS ? GMX_FORCE_NS : 0));
762 /* Clear the unused shake virial and pressure */
763 clear_mat(shake_vir);
764 clear_mat(pres);
766 /* Communicate stuff when parallel */
767 if (PAR(cr) && inputrec->eI != eiNM)
769 wallcycle_start(wcycle, ewcMoveE);
771 global_stat(gstat, cr, enerd, force_vir, shake_vir, mu_tot,
772 inputrec, nullptr, nullptr, nullptr, 1, &terminate,
773 nullptr, FALSE,
774 CGLO_ENERGY |
775 CGLO_PRESSURE |
776 CGLO_CONSTRAINT);
778 wallcycle_stop(wcycle, ewcMoveE);
781 /* Calculate long range corrections to pressure and energy */
782 calc_dispcorr(inputrec, fr, ems->s.box, ems->s.lambda[efptVDW],
783 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
784 enerd->term[F_DISPCORR] = enercorr;
785 enerd->term[F_EPOT] += enercorr;
786 enerd->term[F_PRES] += prescorr;
787 enerd->term[F_DVDL] += dvdlcorr;
789 ems->epot = enerd->term[F_EPOT];
791 if (constr)
793 /* Project out the constraint components of the force */
794 wallcycle_start(wcycle, ewcCONSTR);
795 dvdl_constr = 0;
796 rvec *f_rvec = as_rvec_array(ems->f.data());
797 constrain(nullptr, FALSE, FALSE, constr, &top->idef,
798 inputrec, cr, count, 0, 1.0, mdatoms,
799 as_rvec_array(ems->s.x.data()), f_rvec, f_rvec,
800 fr->bMolPBC, ems->s.box,
801 ems->s.lambda[efptBONDED], &dvdl_constr,
802 nullptr, &shake_vir, nrnb, econqForceDispl);
803 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
804 m_add(force_vir, shake_vir, vir);
805 wallcycle_stop(wcycle, ewcCONSTR);
807 else
809 copy_mat(force_vir, vir);
812 clear_mat(ekin);
813 enerd->term[F_PRES] =
814 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
816 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
818 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
820 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
824 //! Parallel utility summing energies and forces
825 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
826 gmx_mtop_t *top_global,
827 em_state_t *s_min, em_state_t *s_b)
829 t_block *cgs_gl;
830 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
831 double partsum;
832 unsigned char *grpnrFREEZE;
834 if (debug)
836 fprintf(debug, "Doing reorder_partsum\n");
839 const rvec *fm = as_rvec_array(s_min->f.data());
840 const rvec *fb = as_rvec_array(s_b->f.data());
842 cgs_gl = dd_charge_groups_global(cr->dd);
843 index = cgs_gl->index;
845 /* Collect fm in a global vector fmg.
846 * This conflicts with the spirit of domain decomposition,
847 * but to fully optimize this a much more complicated algorithm is required.
849 rvec *fmg;
850 snew(fmg, top_global->natoms);
852 ncg = s_min->s.cg_gl.size();
853 cg_gl = s_min->s.cg_gl.data();
854 i = 0;
855 for (c = 0; c < ncg; c++)
857 cg = cg_gl[c];
858 a0 = index[cg];
859 a1 = index[cg+1];
860 for (a = a0; a < a1; a++)
862 copy_rvec(fm[i], fmg[a]);
863 i++;
866 gmx_sum(top_global->natoms*3, fmg[0], cr);
868 /* Now we will determine the part of the sum for the cgs in state s_b */
869 ncg = s_b->s.cg_gl.size();
870 cg_gl = s_b->s.cg_gl.data();
871 partsum = 0;
872 i = 0;
873 gf = 0;
874 grpnrFREEZE = top_global->groups.grpnr[egcFREEZE];
875 for (c = 0; c < ncg; c++)
877 cg = cg_gl[c];
878 a0 = index[cg];
879 a1 = index[cg+1];
880 for (a = a0; a < a1; a++)
882 if (mdatoms->cFREEZE && grpnrFREEZE)
884 gf = grpnrFREEZE[i];
886 for (m = 0; m < DIM; m++)
888 if (!opts->nFreeze[gf][m])
890 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
893 i++;
897 sfree(fmg);
899 return partsum;
902 //! Print some stuff, like beta, whatever that means.
903 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
904 gmx_mtop_t *top_global,
905 em_state_t *s_min, em_state_t *s_b)
907 double sum;
909 /* This is just the classical Polak-Ribiere calculation of beta;
910 * it looks a bit complicated since we take freeze groups into account,
911 * and might have to sum it in parallel runs.
914 if (!DOMAINDECOMP(cr) ||
915 (s_min->s.ddp_count == cr->dd->ddp_count &&
916 s_b->s.ddp_count == cr->dd->ddp_count))
918 const rvec *fm = as_rvec_array(s_min->f.data());
919 const rvec *fb = as_rvec_array(s_b->f.data());
920 sum = 0;
921 int gf = 0;
922 /* This part of code can be incorrect with DD,
923 * since the atom ordering in s_b and s_min might differ.
925 for (int i = 0; i < mdatoms->homenr; i++)
927 if (mdatoms->cFREEZE)
929 gf = mdatoms->cFREEZE[i];
931 for (int m = 0; m < DIM; m++)
933 if (!opts->nFreeze[gf][m])
935 sum += (fb[i][m] - fm[i][m])*fb[i][m];
940 else
942 /* We need to reorder cgs while summing */
943 sum = reorder_partsum(cr, opts, mdatoms, top_global, s_min, s_b);
945 if (PAR(cr))
947 gmx_sumd(1, &sum, cr);
950 return sum/gmx::square(s_min->fnorm);
953 namespace gmx
956 /*! \brief Do conjugate gradients minimization
957 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
958 int nfile, const t_filenm fnm[],
959 const gmx_output_env_t *oenv, gmx_bool bVerbose,
960 int nstglobalcomm,
961 gmx_vsite_t *vsite, gmx_constr_t constr,
962 int stepout,
963 t_inputrec *inputrec,
964 gmx_mtop_t *top_global, t_fcdata *fcd,
965 t_state *state_global,
966 t_mdatoms *mdatoms,
967 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
968 gmx_edsam_t ed,
969 t_forcerec *fr,
970 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
971 gmx_membed_t gmx_unused *membed,
972 real cpt_period, real max_hours,
973 int imdport,
974 unsigned long Flags,
975 gmx_walltime_accounting_t walltime_accounting)
977 double do_cg(FILE *fplog, t_commrec *cr, const gmx::MDLogger gmx_unused &mdlog,
978 int nfile, const t_filenm fnm[],
979 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
980 int gmx_unused nstglobalcomm,
981 gmx_vsite_t *vsite, gmx_constr_t constr,
982 int gmx_unused stepout,
983 t_inputrec *inputrec,
984 gmx_mtop_t *top_global, t_fcdata *fcd,
985 t_state *state_global,
986 energyhistory_t *energyHistory,
987 t_mdatoms *mdatoms,
988 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
989 gmx_edsam_t gmx_unused ed,
990 t_forcerec *fr,
991 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
992 gmx_membed_t gmx_unused *membed,
993 real gmx_unused cpt_period, real gmx_unused max_hours,
994 int imdport,
995 unsigned long gmx_unused Flags,
996 gmx_walltime_accounting_t walltime_accounting)
998 const char *CG = "Polak-Ribiere Conjugate Gradients";
1000 gmx_localtop_t *top;
1001 gmx_enerdata_t *enerd;
1002 gmx_global_stat_t gstat;
1003 t_graph *graph;
1004 double tmp, minstep;
1005 real stepsize;
1006 real a, b, c, beta = 0.0;
1007 real epot_repl = 0;
1008 real pnorm;
1009 t_mdebin *mdebin;
1010 gmx_bool converged, foundlower;
1011 rvec mu_tot;
1012 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1013 tensor vir, pres;
1014 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1015 gmx_mdoutf_t outf;
1016 int m, step, nminstep;
1018 step = 0;
1020 // Ensure the extra per-atom state array gets allocated
1021 state_global->flags |= (1<<estCGP);
1023 /* Create 4 states on the stack and extract pointers that we will swap */
1024 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1025 em_state_t *s_min = &s0;
1026 em_state_t *s_a = &s1;
1027 em_state_t *s_b = &s2;
1028 em_state_t *s_c = &s3;
1030 /* Init em and store the local state in s_min */
1031 init_em(fplog, CG, cr, inputrec,
1032 state_global, top_global, s_min, &top,
1033 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat,
1034 vsite, constr, nullptr,
1035 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
1037 /* Print to log file */
1038 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1040 /* Max number of steps */
1041 number_steps = inputrec->nsteps;
1043 if (MASTER(cr))
1045 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1047 if (fplog)
1049 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1052 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1053 /* do_force always puts the charge groups in the box and shifts again
1054 * We do not unshift, so molecules are always whole in congrad.c
1056 evaluate_energy(fplog, cr,
1057 top_global, s_min, top,
1058 inputrec, nrnb, wcycle, gstat,
1059 vsite, constr, fcd, graph, mdatoms, fr,
1060 mu_tot, enerd, vir, pres, -1, TRUE);
1061 where();
1063 if (MASTER(cr))
1065 /* Copy stuff to the energy bin for easy printing etc. */
1066 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1067 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1068 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1070 print_ebin_header(fplog, step, step);
1071 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1072 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1074 where();
1076 /* Estimate/guess the initial stepsize */
1077 stepsize = inputrec->em_stepsize/s_min->fnorm;
1079 if (MASTER(cr))
1081 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1082 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1083 s_min->fmax, s_min->a_fmax+1);
1084 fprintf(stderr, " F-Norm = %12.5e\n",
1085 s_min->fnorm/sqrtNumAtoms);
1086 fprintf(stderr, "\n");
1087 /* and copy to the log file too... */
1088 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1089 s_min->fmax, s_min->a_fmax+1);
1090 fprintf(fplog, " F-Norm = %12.5e\n",
1091 s_min->fnorm/sqrtNumAtoms);
1092 fprintf(fplog, "\n");
1094 /* Start the loop over CG steps.
1095 * Each successful step is counted, and we continue until
1096 * we either converge or reach the max number of steps.
1098 converged = FALSE;
1099 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1102 /* start taking steps in a new direction
1103 * First time we enter the routine, beta=0, and the direction is
1104 * simply the negative gradient.
1107 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1108 rvec *pm = as_rvec_array(s_min->s.cg_p.data());
1109 const rvec *sfm = as_rvec_array(s_min->f.data());
1110 double gpa = 0;
1111 int gf = 0;
1112 for (int i = 0; i < mdatoms->homenr; i++)
1114 if (mdatoms->cFREEZE)
1116 gf = mdatoms->cFREEZE[i];
1118 for (m = 0; m < DIM; m++)
1120 if (!inputrec->opts.nFreeze[gf][m])
1122 pm[i][m] = sfm[i][m] + beta*pm[i][m];
1123 gpa -= pm[i][m]*sfm[i][m];
1124 /* f is negative gradient, thus the sign */
1126 else
1128 pm[i][m] = 0;
1133 /* Sum the gradient along the line across CPUs */
1134 if (PAR(cr))
1136 gmx_sumd(1, &gpa, cr);
1139 /* Calculate the norm of the search vector */
1140 get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1142 /* Just in case stepsize reaches zero due to numerical precision... */
1143 if (stepsize <= 0)
1145 stepsize = inputrec->em_stepsize/pnorm;
1149 * Double check the value of the derivative in the search direction.
1150 * If it is positive it must be due to the old information in the
1151 * CG formula, so just remove that and start over with beta=0.
1152 * This corresponds to a steepest descent step.
1154 if (gpa > 0)
1156 beta = 0;
1157 step--; /* Don't count this step since we are restarting */
1158 continue; /* Go back to the beginning of the big for-loop */
1161 /* Calculate minimum allowed stepsize, before the average (norm)
1162 * relative change in coordinate is smaller than precision
1164 minstep = 0;
1165 for (int i = 0; i < mdatoms->homenr; i++)
1167 for (m = 0; m < DIM; m++)
1169 tmp = fabs(s_min->s.x[i][m]);
1170 if (tmp < 1.0)
1172 tmp = 1.0;
1174 tmp = pm[i][m]/tmp;
1175 minstep += tmp*tmp;
1178 /* Add up from all CPUs */
1179 if (PAR(cr))
1181 gmx_sumd(1, &minstep, cr);
1184 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1186 if (stepsize < minstep)
1188 converged = TRUE;
1189 break;
1192 /* Write coordinates if necessary */
1193 do_x = do_per_step(step, inputrec->nstxout);
1194 do_f = do_per_step(step, inputrec->nstfout);
1196 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
1197 top_global, inputrec, step,
1198 s_min, state_global, energyHistory);
1200 /* Take a step downhill.
1201 * In theory, we should minimize the function along this direction.
1202 * That is quite possible, but it turns out to take 5-10 function evaluations
1203 * for each line. However, we dont really need to find the exact minimum -
1204 * it is much better to start a new CG step in a modified direction as soon
1205 * as we are close to it. This will save a lot of energy evaluations.
1207 * In practice, we just try to take a single step.
1208 * If it worked (i.e. lowered the energy), we increase the stepsize but
1209 * the continue straight to the next CG step without trying to find any minimum.
1210 * If it didn't work (higher energy), there must be a minimum somewhere between
1211 * the old position and the new one.
1213 * Due to the finite numerical accuracy, it turns out that it is a good idea
1214 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1215 * This leads to lower final energies in the tests I've done. / Erik
1217 s_a->epot = s_min->epot;
1218 a = 0.0;
1219 c = a + stepsize; /* reference position along line is zero */
1221 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1223 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1224 s_min, top, mdatoms, fr, vsite, constr,
1225 nrnb, wcycle);
1228 /* Take a trial step (new coords in s_c) */
1229 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, &s_min->s.cg_p, s_c,
1230 constr, top, nrnb, wcycle, -1);
1232 neval++;
1233 /* Calculate energy for the trial step */
1234 evaluate_energy(fplog, cr,
1235 top_global, s_c, top,
1236 inputrec, nrnb, wcycle, gstat,
1237 vsite, constr, fcd, graph, mdatoms, fr,
1238 mu_tot, enerd, vir, pres, -1, FALSE);
1240 /* Calc derivative along line */
1241 const rvec *pc = as_rvec_array(s_c->s.cg_p.data());
1242 const rvec *sfc = as_rvec_array(s_c->f.data());
1243 double gpc = 0;
1244 for (int i = 0; i < mdatoms->homenr; i++)
1246 for (m = 0; m < DIM; m++)
1248 gpc -= pc[i][m]*sfc[i][m]; /* f is negative gradient, thus the sign */
1251 /* Sum the gradient along the line across CPUs */
1252 if (PAR(cr))
1254 gmx_sumd(1, &gpc, cr);
1257 /* This is the max amount of increase in energy we tolerate */
1258 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1260 /* Accept the step if the energy is lower, or if it is not significantly higher
1261 * and the line derivative is still negative.
1263 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1265 foundlower = TRUE;
1266 /* Great, we found a better energy. Increase step for next iteration
1267 * if we are still going down, decrease it otherwise
1269 if (gpc < 0)
1271 stepsize *= 1.618034; /* The golden section */
1273 else
1275 stepsize *= 0.618034; /* 1/golden section */
1278 else
1280 /* New energy is the same or higher. We will have to do some work
1281 * to find a smaller value in the interval. Take smaller step next time!
1283 foundlower = FALSE;
1284 stepsize *= 0.618034;
1290 /* OK, if we didn't find a lower value we will have to locate one now - there must
1291 * be one in the interval [a=0,c].
1292 * The same thing is valid here, though: Don't spend dozens of iterations to find
1293 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1294 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1296 * I also have a safeguard for potentially really pathological functions so we never
1297 * take more than 20 steps before we give up ...
1299 * If we already found a lower value we just skip this step and continue to the update.
1301 double gpb;
1302 if (!foundlower)
1304 nminstep = 0;
1308 /* Select a new trial point.
1309 * If the derivatives at points a & c have different sign we interpolate to zero,
1310 * otherwise just do a bisection.
1312 if (gpa < 0 && gpc > 0)
1314 b = a + gpa*(a-c)/(gpc-gpa);
1316 else
1318 b = 0.5*(a+c);
1321 /* safeguard if interpolation close to machine accuracy causes errors:
1322 * never go outside the interval
1324 if (b <= a || b >= c)
1326 b = 0.5*(a+c);
1329 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1331 /* Reload the old state */
1332 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1333 s_min, top, mdatoms, fr, vsite, constr,
1334 nrnb, wcycle);
1337 /* Take a trial step to this new point - new coords in s_b */
1338 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, &s_min->s.cg_p, s_b,
1339 constr, top, nrnb, wcycle, -1);
1341 neval++;
1342 /* Calculate energy for the trial step */
1343 evaluate_energy(fplog, cr,
1344 top_global, s_b, top,
1345 inputrec, nrnb, wcycle, gstat,
1346 vsite, constr, fcd, graph, mdatoms, fr,
1347 mu_tot, enerd, vir, pres, -1, FALSE);
1349 /* p does not change within a step, but since the domain decomposition
1350 * might change, we have to use cg_p of s_b here.
1352 const rvec *pb = as_rvec_array(s_b->s.cg_p.data());
1353 const rvec *sfb = as_rvec_array(s_b->f.data());
1354 gpb = 0;
1355 for (int i = 0; i < mdatoms->homenr; i++)
1357 for (m = 0; m < DIM; m++)
1359 gpb -= pb[i][m]*sfb[i][m]; /* f is negative gradient, thus the sign */
1362 /* Sum the gradient along the line across CPUs */
1363 if (PAR(cr))
1365 gmx_sumd(1, &gpb, cr);
1368 if (debug)
1370 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1371 s_a->epot, s_b->epot, s_c->epot, gpb);
1374 epot_repl = s_b->epot;
1376 /* Keep one of the intervals based on the value of the derivative at the new point */
1377 if (gpb > 0)
1379 /* Replace c endpoint with b */
1380 swap_em_state(&s_b, &s_c);
1381 c = b;
1382 gpc = gpb;
1384 else
1386 /* Replace a endpoint with b */
1387 swap_em_state(&s_b, &s_a);
1388 a = b;
1389 gpa = gpb;
1393 * Stop search as soon as we find a value smaller than the endpoints.
1394 * Never run more than 20 steps, no matter what.
1396 nminstep++;
1398 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1399 (nminstep < 20));
1401 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1402 nminstep >= 20)
1404 /* OK. We couldn't find a significantly lower energy.
1405 * If beta==0 this was steepest descent, and then we give up.
1406 * If not, set beta=0 and restart with steepest descent before quitting.
1408 if (beta == 0.0)
1410 /* Converged */
1411 converged = TRUE;
1412 break;
1414 else
1416 /* Reset memory before giving up */
1417 beta = 0.0;
1418 continue;
1422 /* Select min energy state of A & C, put the best in B.
1424 if (s_c->epot < s_a->epot)
1426 if (debug)
1428 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1429 s_c->epot, s_a->epot);
1431 swap_em_state(&s_b, &s_c);
1432 gpb = gpc;
1434 else
1436 if (debug)
1438 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1439 s_a->epot, s_c->epot);
1441 swap_em_state(&s_b, &s_a);
1442 gpb = gpa;
1446 else
1448 if (debug)
1450 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1451 s_c->epot);
1453 swap_em_state(&s_b, &s_c);
1454 gpb = gpc;
1457 /* new search direction */
1458 /* beta = 0 means forget all memory and restart with steepest descents. */
1459 if (nstcg && ((step % nstcg) == 0))
1461 beta = 0.0;
1463 else
1465 /* s_min->fnorm cannot be zero, because then we would have converged
1466 * and broken out.
1469 /* Polak-Ribiere update.
1470 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1472 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1474 /* Limit beta to prevent oscillations */
1475 if (fabs(beta) > 5.0)
1477 beta = 0.0;
1481 /* update positions */
1482 swap_em_state(&s_min, &s_b);
1483 gpa = gpb;
1485 /* Print it if necessary */
1486 if (MASTER(cr))
1488 if (bVerbose)
1490 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1491 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1492 step, s_min->epot, s_min->fnorm/sqrtNumAtoms,
1493 s_min->fmax, s_min->a_fmax+1);
1494 fflush(stderr);
1496 /* Store the new (lower) energies */
1497 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1498 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1499 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1501 do_log = do_per_step(step, inputrec->nstlog);
1502 do_ene = do_per_step(step, inputrec->nstenergy);
1504 /* Prepare IMD energy record, if bIMD is TRUE. */
1505 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1507 if (do_log)
1509 print_ebin_header(fplog, step, step);
1511 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1512 do_log ? fplog : nullptr, step, step, eprNORMAL,
1513 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1516 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1517 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, as_rvec_array(state_global->x.data()), inputrec, 0, wcycle) && MASTER(cr))
1519 IMD_send_positions(inputrec->imd);
1522 /* Stop when the maximum force lies below tolerance.
1523 * If we have reached machine precision, converged is already set to true.
1525 converged = converged || (s_min->fmax < inputrec->em_tol);
1527 } /* End of the loop */
1529 /* IMD cleanup, if bIMD is TRUE. */
1530 IMD_finalize(inputrec->bIMD, inputrec->imd);
1532 if (converged)
1534 step--; /* we never took that last step in this case */
1537 if (s_min->fmax > inputrec->em_tol)
1539 if (MASTER(cr))
1541 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1542 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1544 converged = FALSE;
1547 if (MASTER(cr))
1549 /* If we printed energy and/or logfile last step (which was the last step)
1550 * we don't have to do it again, but otherwise print the final values.
1552 if (!do_log)
1554 /* Write final value to log since we didn't do anything the last step */
1555 print_ebin_header(fplog, step, step);
1557 if (!do_ene || !do_log)
1559 /* Write final energy file entries */
1560 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1561 !do_log ? fplog : nullptr, step, step, eprNORMAL,
1562 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1566 /* Print some stuff... */
1567 if (MASTER(cr))
1569 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1572 /* IMPORTANT!
1573 * For accurate normal mode calculation it is imperative that we
1574 * store the last conformation into the full precision binary trajectory.
1576 * However, we should only do it if we did NOT already write this step
1577 * above (which we did if do_x or do_f was true).
1579 do_x = !do_per_step(step, inputrec->nstxout);
1580 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1582 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1583 top_global, inputrec, step,
1584 s_min, state_global, energyHistory);
1587 if (MASTER(cr))
1589 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1590 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1591 s_min, sqrtNumAtoms);
1592 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1593 s_min, sqrtNumAtoms);
1595 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1598 finish_em(cr, outf, inputrec, walltime_accounting, wcycle);
1600 /* To print the actual number of steps we needed somewhere */
1601 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1603 return 0;
1604 } /* That's all folks */
1607 /*! \brief Do L-BFGS conjugate gradients minimization
1608 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
1609 int nfile, const t_filenm fnm[],
1610 const gmx_output_env_t *oenv, gmx_bool bVerbose,
1611 int nstglobalcomm,
1612 gmx_vsite_t *vsite, gmx_constr_t constr,
1613 int stepout,
1614 t_inputrec *inputrec,
1615 gmx_mtop_t *top_global, t_fcdata *fcd,
1616 t_state *state_global,
1617 t_mdatoms *mdatoms,
1618 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1619 gmx_edsam_t ed,
1620 t_forcerec *fr,
1621 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
1622 real cpt_period, real max_hours,
1623 int imdport,
1624 unsigned long Flags,
1625 gmx_walltime_accounting_t walltime_accounting)
1627 double do_lbfgs(FILE *fplog, t_commrec *cr, const gmx::MDLogger gmx_unused &mdlog,
1628 int nfile, const t_filenm fnm[],
1629 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
1630 int gmx_unused nstglobalcomm,
1631 gmx_vsite_t *vsite, gmx_constr_t constr,
1632 int gmx_unused stepout,
1633 t_inputrec *inputrec,
1634 gmx_mtop_t *top_global, t_fcdata *fcd,
1635 t_state *state_global,
1636 energyhistory_t *energyHistory,
1637 t_mdatoms *mdatoms,
1638 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1639 gmx_edsam_t gmx_unused ed,
1640 t_forcerec *fr,
1641 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1642 gmx_membed_t gmx_unused *membed,
1643 real gmx_unused cpt_period, real gmx_unused max_hours,
1644 int imdport,
1645 unsigned long gmx_unused Flags,
1646 gmx_walltime_accounting_t walltime_accounting)
1648 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1649 em_state_t ems;
1650 gmx_localtop_t *top;
1651 gmx_enerdata_t *enerd;
1652 gmx_global_stat_t gstat;
1653 t_graph *graph;
1654 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1655 double stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1656 real *rho, *alpha, *p, *s, **dx, **dg;
1657 real a, b, c, maxdelta, delta;
1658 real diag, Epot0;
1659 real dgdx, dgdg, sq, yr, beta;
1660 t_mdebin *mdebin;
1661 gmx_bool converged;
1662 rvec mu_tot;
1663 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1664 tensor vir, pres;
1665 int start, end, number_steps;
1666 gmx_mdoutf_t outf;
1667 int i, k, m, n, gf, step;
1668 int mdof_flags;
1670 if (PAR(cr))
1672 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1675 if (nullptr != constr)
1677 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).");
1680 n = 3*state_global->natoms;
1681 nmaxcorr = inputrec->nbfgscorr;
1683 snew(frozen, n);
1685 snew(p, n);
1686 snew(rho, nmaxcorr);
1687 snew(alpha, nmaxcorr);
1689 snew(dx, nmaxcorr);
1690 for (i = 0; i < nmaxcorr; i++)
1692 snew(dx[i], n);
1695 snew(dg, nmaxcorr);
1696 for (i = 0; i < nmaxcorr; i++)
1698 snew(dg[i], n);
1701 step = 0;
1702 neval = 0;
1704 /* Init em */
1705 init_em(fplog, LBFGS, cr, inputrec,
1706 state_global, top_global, &ems, &top,
1707 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat,
1708 vsite, constr, nullptr,
1709 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
1711 start = 0;
1712 end = mdatoms->homenr;
1714 /* We need 4 working states */
1715 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1716 em_state_t *sa = &s0;
1717 em_state_t *sb = &s1;
1718 em_state_t *sc = &s2;
1719 em_state_t *last = &s3;
1720 /* Initialize by copying the state from ems (we could skip x and f here) */
1721 *sa = ems;
1722 *sb = ems;
1723 *sc = ems;
1725 /* Print to log file */
1726 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1728 do_log = do_ene = do_x = do_f = TRUE;
1730 /* Max number of steps */
1731 number_steps = inputrec->nsteps;
1733 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1734 gf = 0;
1735 for (i = start; i < end; i++)
1737 if (mdatoms->cFREEZE)
1739 gf = mdatoms->cFREEZE[i];
1741 for (m = 0; m < DIM; m++)
1743 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1746 if (MASTER(cr))
1748 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1750 if (fplog)
1752 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1755 if (vsite)
1757 construct_vsites(vsite, as_rvec_array(state_global->x.data()), 1, nullptr,
1758 top->idef.iparams, top->idef.il,
1759 fr->ePBC, fr->bMolPBC, cr, state_global->box);
1762 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1763 /* do_force always puts the charge groups in the box and shifts again
1764 * We do not unshift, so molecules are always whole
1766 neval++;
1767 evaluate_energy(fplog, cr,
1768 top_global, &ems, top,
1769 inputrec, nrnb, wcycle, gstat,
1770 vsite, constr, fcd, graph, mdatoms, fr,
1771 mu_tot, enerd, vir, pres, -1, TRUE);
1772 where();
1774 if (MASTER(cr))
1776 /* Copy stuff to the energy bin for easy printing etc. */
1777 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1778 mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
1779 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1781 print_ebin_header(fplog, step, step);
1782 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1783 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1785 where();
1787 /* Set the initial step.
1788 * since it will be multiplied by the non-normalized search direction
1789 * vector (force vector the first time), we scale it by the
1790 * norm of the force.
1793 if (MASTER(cr))
1795 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1796 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1797 fprintf(stderr, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1798 fprintf(stderr, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1799 fprintf(stderr, "\n");
1800 /* and copy to the log file too... */
1801 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1802 fprintf(fplog, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1803 fprintf(fplog, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1804 fprintf(fplog, "\n");
1807 // Point is an index to the memory of search directions, where 0 is the first one.
1808 point = 0;
1810 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1811 real *fInit = static_cast<real *>(as_rvec_array(ems.f.data())[0]);
1812 for (i = 0; i < n; i++)
1814 if (!frozen[i])
1816 dx[point][i] = fInit[i]; /* Initial search direction */
1818 else
1820 dx[point][i] = 0;
1824 // Stepsize will be modified during the search, and actually it is not critical
1825 // (the main efficiency in the algorithm comes from changing directions), but
1826 // we still need an initial value, so estimate it as the inverse of the norm
1827 // so we take small steps where the potential fluctuates a lot.
1828 stepsize = 1.0/ems.fnorm;
1830 /* Start the loop over BFGS steps.
1831 * Each successful step is counted, and we continue until
1832 * we either converge or reach the max number of steps.
1835 ncorr = 0;
1837 /* Set the gradient from the force */
1838 converged = FALSE;
1839 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1842 /* Write coordinates if necessary */
1843 do_x = do_per_step(step, inputrec->nstxout);
1844 do_f = do_per_step(step, inputrec->nstfout);
1846 mdof_flags = 0;
1847 if (do_x)
1849 mdof_flags |= MDOF_X;
1852 if (do_f)
1854 mdof_flags |= MDOF_F;
1857 if (inputrec->bIMD)
1859 mdof_flags |= MDOF_IMD;
1862 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1863 top_global, step, (real)step, &ems.s, state_global, energyHistory, &ems.f);
1865 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1867 /* make s a pointer to current search direction - point=0 first time we get here */
1868 s = dx[point];
1870 real *xx = static_cast<real *>(as_rvec_array(ems.s.x.data())[0]);
1871 real *ff = static_cast<real *>(as_rvec_array(ems.f.data())[0]);
1873 // calculate line gradient in position A
1874 for (gpa = 0, i = 0; i < n; i++)
1876 gpa -= s[i]*ff[i];
1879 /* Calculate minimum allowed stepsize along the line, before the average (norm)
1880 * relative change in coordinate is smaller than precision
1882 for (minstep = 0, i = 0; i < n; i++)
1884 tmp = fabs(xx[i]);
1885 if (tmp < 1.0)
1887 tmp = 1.0;
1889 tmp = s[i]/tmp;
1890 minstep += tmp*tmp;
1892 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1894 if (stepsize < minstep)
1896 converged = TRUE;
1897 break;
1900 // Before taking any steps along the line, store the old position
1901 *last = ems;
1902 real *lastx = static_cast<real *>(as_rvec_array(last->s.x.data())[0]);
1903 real *lastf = static_cast<real *>(as_rvec_array(last->f.data())[0]);
1904 Epot0 = ems.epot;
1906 *sa = ems;
1908 /* Take a step downhill.
1909 * In theory, we should find the actual minimum of the function in this
1910 * direction, somewhere along the line.
1911 * That is quite possible, but it turns out to take 5-10 function evaluations
1912 * for each line. However, we dont really need to find the exact minimum -
1913 * it is much better to start a new BFGS step in a modified direction as soon
1914 * as we are close to it. This will save a lot of energy evaluations.
1916 * In practice, we just try to take a single step.
1917 * If it worked (i.e. lowered the energy), we increase the stepsize but
1918 * continue straight to the next BFGS step without trying to find any minimum,
1919 * i.e. we change the search direction too. If the line was smooth, it is
1920 * likely we are in a smooth region, and then it makes sense to take longer
1921 * steps in the modified search direction too.
1923 * If it didn't work (higher energy), there must be a minimum somewhere between
1924 * the old position and the new one. Then we need to start by finding a lower
1925 * value before we change search direction. Since the energy was apparently
1926 * quite rough, we need to decrease the step size.
1928 * Due to the finite numerical accuracy, it turns out that it is a good idea
1929 * to accept a SMALL increase in energy, if the derivative is still downhill.
1930 * This leads to lower final energies in the tests I've done. / Erik
1933 // State "A" is the first position along the line.
1934 // reference position along line is initially zero
1935 a = 0.0;
1937 // Check stepsize first. We do not allow displacements
1938 // larger than emstep.
1942 // Pick a new position C by adding stepsize to A.
1943 c = a + stepsize;
1945 // Calculate what the largest change in any individual coordinate
1946 // would be (translation along line * gradient along line)
1947 maxdelta = 0;
1948 for (i = 0; i < n; i++)
1950 delta = c*s[i];
1951 if (delta > maxdelta)
1953 maxdelta = delta;
1956 // If any displacement is larger than the stepsize limit, reduce the step
1957 if (maxdelta > inputrec->em_stepsize)
1959 stepsize *= 0.1;
1962 while (maxdelta > inputrec->em_stepsize);
1964 // Take a trial step and move the coordinate array xc[] to position C
1965 real *xc = static_cast<real *>(as_rvec_array(sc->s.x.data())[0]);
1966 for (i = 0; i < n; i++)
1968 xc[i] = lastx[i] + c*s[i];
1971 neval++;
1972 // Calculate energy for the trial step in position C
1973 evaluate_energy(fplog, cr,
1974 top_global, sc, top,
1975 inputrec, nrnb, wcycle, gstat,
1976 vsite, constr, fcd, graph, mdatoms, fr,
1977 mu_tot, enerd, vir, pres, step, FALSE);
1979 // Calc line gradient in position C
1980 real *fc = static_cast<real *>(as_rvec_array(sc->f.data())[0]);
1981 for (gpc = 0, i = 0; i < n; i++)
1983 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1985 /* Sum the gradient along the line across CPUs */
1986 if (PAR(cr))
1988 gmx_sumd(1, &gpc, cr);
1991 // This is the max amount of increase in energy we tolerate.
1992 // By allowing VERY small changes (close to numerical precision) we
1993 // frequently find even better (lower) final energies.
1994 tmp = sqrt(GMX_REAL_EPS)*fabs(sa->epot);
1996 // Accept the step if the energy is lower in the new position C (compared to A),
1997 // or if it is not significantly higher and the line derivative is still negative.
1998 if (sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp)))
2000 // Great, we found a better energy. We no longer try to alter the
2001 // stepsize, but simply accept this new better position. The we select a new
2002 // search direction instead, which will be much more efficient than continuing
2003 // to take smaller steps along a line. Set fnorm based on the new C position,
2004 // which will be used to update the stepsize to 1/fnorm further down.
2005 foundlower = TRUE;
2007 else
2009 // If we got here, the energy is NOT lower in point C, i.e. it will be the same
2010 // or higher than in point A. In this case it is pointless to move to point C,
2011 // so we will have to do more iterations along the same line to find a smaller
2012 // value in the interval [A=0.0,C].
2013 // Here, A is still 0.0, but that will change when we do a search in the interval
2014 // [0.0,C] below. That search we will do by interpolation or bisection rather
2015 // than with the stepsize, so no need to modify it. For the next search direction
2016 // it will be reset to 1/fnorm anyway.
2017 foundlower = FALSE;
2020 if (!foundlower)
2022 // OK, if we didn't find a lower value we will have to locate one now - there must
2023 // be one in the interval [a,c].
2024 // The same thing is valid here, though: Don't spend dozens of iterations to find
2025 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2026 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2027 // I also have a safeguard for potentially really pathological functions so we never
2028 // take more than 20 steps before we give up.
2029 // If we already found a lower value we just skip this step and continue to the update.
2030 real fnorm = 0;
2031 nminstep = 0;
2034 // Select a new trial point B in the interval [A,C].
2035 // If the derivatives at points a & c have different sign we interpolate to zero,
2036 // otherwise just do a bisection since there might be multiple minima/maxima
2037 // inside the interval.
2038 if (gpa < 0 && gpc > 0)
2040 b = a + gpa*(a-c)/(gpc-gpa);
2042 else
2044 b = 0.5*(a+c);
2047 /* safeguard if interpolation close to machine accuracy causes errors:
2048 * never go outside the interval
2050 if (b <= a || b >= c)
2052 b = 0.5*(a+c);
2055 // Take a trial step to point B
2056 real *xb = static_cast<real *>(as_rvec_array(sb->s.x.data())[0]);
2057 for (i = 0; i < n; i++)
2059 xb[i] = lastx[i] + b*s[i];
2062 neval++;
2063 // Calculate energy for the trial step in point B
2064 evaluate_energy(fplog, cr,
2065 top_global, sb, top,
2066 inputrec, nrnb, wcycle, gstat,
2067 vsite, constr, fcd, graph, mdatoms, fr,
2068 mu_tot, enerd, vir, pres, step, FALSE);
2069 fnorm = sb->fnorm;
2071 // Calculate gradient in point B
2072 real *fb = static_cast<real *>(as_rvec_array(sb->f.data())[0]);
2073 for (gpb = 0, i = 0; i < n; i++)
2075 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2078 /* Sum the gradient along the line across CPUs */
2079 if (PAR(cr))
2081 gmx_sumd(1, &gpb, cr);
2084 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2085 // at the new point B, and rename the endpoints of this new interval A and C.
2086 if (gpb > 0)
2088 /* Replace c endpoint with b */
2089 c = b;
2090 /* swap states b and c */
2091 swap_em_state(&sb, &sc);
2093 else
2095 /* Replace a endpoint with b */
2096 a = b;
2097 /* swap states a and b */
2098 swap_em_state(&sa, &sb);
2102 * Stop search as soon as we find a value smaller than the endpoints,
2103 * or if the tolerance is below machine precision.
2104 * Never run more than 20 steps, no matter what.
2106 nminstep++;
2108 while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2110 if (fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2112 /* OK. We couldn't find a significantly lower energy.
2113 * If ncorr==0 this was steepest descent, and then we give up.
2114 * If not, reset memory to restart as steepest descent before quitting.
2116 if (ncorr == 0)
2118 /* Converged */
2119 converged = TRUE;
2120 break;
2122 else
2124 /* Reset memory */
2125 ncorr = 0;
2126 /* Search in gradient direction */
2127 for (i = 0; i < n; i++)
2129 dx[point][i] = ff[i];
2131 /* Reset stepsize */
2132 stepsize = 1.0/fnorm;
2133 continue;
2137 /* Select min energy state of A & C, put the best in xx/ff/Epot
2139 if (sc->epot < sa->epot)
2141 /* Use state C */
2142 ems = *sc;
2143 step_taken = c;
2145 else
2147 /* Use state A */
2148 ems = *sa;
2149 step_taken = a;
2153 else
2155 /* found lower */
2156 /* Use state C */
2157 ems = *sc;
2158 step_taken = c;
2161 /* Update the memory information, and calculate a new
2162 * approximation of the inverse hessian
2165 /* Have new data in Epot, xx, ff */
2166 if (ncorr < nmaxcorr)
2168 ncorr++;
2171 for (i = 0; i < n; i++)
2173 dg[point][i] = lastf[i]-ff[i];
2174 dx[point][i] *= step_taken;
2177 dgdg = 0;
2178 dgdx = 0;
2179 for (i = 0; i < n; i++)
2181 dgdg += dg[point][i]*dg[point][i];
2182 dgdx += dg[point][i]*dx[point][i];
2185 diag = dgdx/dgdg;
2187 rho[point] = 1.0/dgdx;
2188 point++;
2190 if (point >= nmaxcorr)
2192 point = 0;
2195 /* Update */
2196 for (i = 0; i < n; i++)
2198 p[i] = ff[i];
2201 cp = point;
2203 /* Recursive update. First go back over the memory points */
2204 for (k = 0; k < ncorr; k++)
2206 cp--;
2207 if (cp < 0)
2209 cp = ncorr-1;
2212 sq = 0;
2213 for (i = 0; i < n; i++)
2215 sq += dx[cp][i]*p[i];
2218 alpha[cp] = rho[cp]*sq;
2220 for (i = 0; i < n; i++)
2222 p[i] -= alpha[cp]*dg[cp][i];
2226 for (i = 0; i < n; i++)
2228 p[i] *= diag;
2231 /* And then go forward again */
2232 for (k = 0; k < ncorr; k++)
2234 yr = 0;
2235 for (i = 0; i < n; i++)
2237 yr += p[i]*dg[cp][i];
2240 beta = rho[cp]*yr;
2241 beta = alpha[cp]-beta;
2243 for (i = 0; i < n; i++)
2245 p[i] += beta*dx[cp][i];
2248 cp++;
2249 if (cp >= ncorr)
2251 cp = 0;
2255 for (i = 0; i < n; i++)
2257 if (!frozen[i])
2259 dx[point][i] = p[i];
2261 else
2263 dx[point][i] = 0;
2267 /* Print it if necessary */
2268 if (MASTER(cr))
2270 if (bVerbose)
2272 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2273 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2274 step, ems.epot, ems.fnorm/sqrtNumAtoms, ems.fmax, ems.a_fmax + 1);
2275 fflush(stderr);
2277 /* Store the new (lower) energies */
2278 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2279 mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
2280 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2281 do_log = do_per_step(step, inputrec->nstlog);
2282 do_ene = do_per_step(step, inputrec->nstenergy);
2283 if (do_log)
2285 print_ebin_header(fplog, step, step);
2287 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2288 do_log ? fplog : nullptr, step, step, eprNORMAL,
2289 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2292 /* Send x and E to IMD client, if bIMD is TRUE. */
2293 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, as_rvec_array(state_global->x.data()), inputrec, 0, wcycle) && MASTER(cr))
2295 IMD_send_positions(inputrec->imd);
2298 // Reset stepsize in we are doing more iterations
2299 stepsize = 1.0/ems.fnorm;
2301 /* Stop when the maximum force lies below tolerance.
2302 * If we have reached machine precision, converged is already set to true.
2304 converged = converged || (ems.fmax < inputrec->em_tol);
2306 } /* End of the loop */
2308 /* IMD cleanup, if bIMD is TRUE. */
2309 IMD_finalize(inputrec->bIMD, inputrec->imd);
2311 if (converged)
2313 step--; /* we never took that last step in this case */
2316 if (ems.fmax > inputrec->em_tol)
2318 if (MASTER(cr))
2320 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2321 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2323 converged = FALSE;
2326 /* If we printed energy and/or logfile last step (which was the last step)
2327 * we don't have to do it again, but otherwise print the final values.
2329 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2331 print_ebin_header(fplog, step, step);
2333 if (!do_ene || !do_log) /* Write final energy file entries */
2335 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2336 !do_log ? fplog : nullptr, step, step, eprNORMAL,
2337 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2340 /* Print some stuff... */
2341 if (MASTER(cr))
2343 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2346 /* IMPORTANT!
2347 * For accurate normal mode calculation it is imperative that we
2348 * store the last conformation into the full precision binary trajectory.
2350 * However, we should only do it if we did NOT already write this step
2351 * above (which we did if do_x or do_f was true).
2353 do_x = !do_per_step(step, inputrec->nstxout);
2354 do_f = !do_per_step(step, inputrec->nstfout);
2355 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2356 top_global, inputrec, step,
2357 &ems, state_global, energyHistory);
2359 if (MASTER(cr))
2361 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2362 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2363 number_steps, &ems, sqrtNumAtoms);
2364 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2365 number_steps, &ems, sqrtNumAtoms);
2367 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2370 finish_em(cr, outf, inputrec, walltime_accounting, wcycle);
2372 /* To print the actual number of steps we needed somewhere */
2373 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2375 return 0;
2376 } /* That's all folks */
2378 /*! \brief Do steepest descents minimization
2379 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
2380 int nfile, const t_filenm fnm[],
2381 const gmx_output_env_t *oenv, gmx_bool bVerbose,
2382 int nstglobalcomm,
2383 gmx_vsite_t *vsite, gmx_constr_t constr,
2384 int stepout,
2385 t_inputrec *inputrec,
2386 gmx_mtop_t *top_global, t_fcdata *fcd,
2387 t_state *state_global,
2388 t_mdatoms *mdatoms,
2389 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2390 gmx_edsam_t ed,
2391 t_forcerec *fr,
2392 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
2393 real cpt_period, real max_hours,
2394 int imdport,
2395 unsigned long Flags,
2396 gmx_walltime_accounting_t walltime_accounting)
2398 double do_steep(FILE *fplog, t_commrec *cr, const gmx::MDLogger gmx_unused &mdlog,
2399 int nfile, const t_filenm fnm[],
2400 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
2401 int gmx_unused nstglobalcomm,
2402 gmx_vsite_t *vsite, gmx_constr_t constr,
2403 int gmx_unused stepout,
2404 t_inputrec *inputrec,
2405 gmx_mtop_t *top_global, t_fcdata *fcd,
2406 t_state *state_global,
2407 energyhistory_t *energyHistory,
2408 t_mdatoms *mdatoms,
2409 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2410 gmx_edsam_t gmx_unused ed,
2411 t_forcerec *fr,
2412 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2413 gmx_membed_t gmx_unused *membed,
2414 real gmx_unused cpt_period, real gmx_unused max_hours,
2415 int imdport,
2416 unsigned long gmx_unused Flags,
2417 gmx_walltime_accounting_t walltime_accounting)
2419 const char *SD = "Steepest Descents";
2420 gmx_localtop_t *top;
2421 gmx_enerdata_t *enerd;
2422 gmx_global_stat_t gstat;
2423 t_graph *graph;
2424 real stepsize;
2425 real ustep;
2426 gmx_mdoutf_t outf;
2427 t_mdebin *mdebin;
2428 gmx_bool bDone, bAbort, do_x, do_f;
2429 tensor vir, pres;
2430 rvec mu_tot;
2431 int nsteps;
2432 int count = 0;
2433 int steps_accepted = 0;
2435 /* Create 2 states on the stack and extract pointers that we will swap */
2436 em_state_t s0 {}, s1 {};
2437 em_state_t *s_min = &s0;
2438 em_state_t *s_try = &s1;
2440 /* Init em and store the local state in s_try */
2441 init_em(fplog, SD, cr, inputrec,
2442 state_global, top_global, s_try, &top,
2443 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat,
2444 vsite, constr, nullptr,
2445 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
2447 /* Print to log file */
2448 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2450 /* Set variables for stepsize (in nm). This is the largest
2451 * step that we are going to make in any direction.
2453 ustep = inputrec->em_stepsize;
2454 stepsize = 0;
2456 /* Max number of steps */
2457 nsteps = inputrec->nsteps;
2459 if (MASTER(cr))
2461 /* Print to the screen */
2462 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2464 if (fplog)
2466 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2469 /**** HERE STARTS THE LOOP ****
2470 * count is the counter for the number of steps
2471 * bDone will be TRUE when the minimization has converged
2472 * bAbort will be TRUE when nsteps steps have been performed or when
2473 * the stepsize becomes smaller than is reasonable for machine precision
2475 count = 0;
2476 bDone = FALSE;
2477 bAbort = FALSE;
2478 while (!bDone && !bAbort)
2480 bAbort = (nsteps >= 0) && (count == nsteps);
2482 /* set new coordinates, except for first step */
2483 bool validStep = true;
2484 if (count > 0)
2486 validStep =
2487 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2488 s_min, stepsize, &s_min->f, s_try,
2489 constr, top, nrnb, wcycle, count);
2492 if (validStep)
2494 evaluate_energy(fplog, cr,
2495 top_global, s_try, top,
2496 inputrec, nrnb, wcycle, gstat,
2497 vsite, constr, fcd, graph, mdatoms, fr,
2498 mu_tot, enerd, vir, pres, count, count == 0);
2500 else
2502 // Signal constraint error during stepping with energy=inf
2503 s_try->epot = std::numeric_limits<real>::infinity();
2506 if (MASTER(cr))
2508 print_ebin_header(fplog, count, count);
2511 if (count == 0)
2513 s_min->epot = s_try->epot;
2516 /* Print it if necessary */
2517 if (MASTER(cr))
2519 if (bVerbose)
2521 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2522 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2523 ( (count == 0) || (s_try->epot < s_min->epot) ) ? '\n' : '\r');
2524 fflush(stderr);
2527 if ( (count == 0) || (s_try->epot < s_min->epot) )
2529 /* Store the new (lower) energies */
2530 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2531 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2532 s_try->s.box, nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2534 /* Prepare IMD energy record, if bIMD is TRUE. */
2535 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2537 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2538 do_per_step(steps_accepted, inputrec->nstdisreout),
2539 do_per_step(steps_accepted, inputrec->nstorireout),
2540 fplog, count, count, eprNORMAL,
2541 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2542 fflush(fplog);
2546 /* Now if the new energy is smaller than the previous...
2547 * or if this is the first step!
2548 * or if we did random steps!
2551 if ( (count == 0) || (s_try->epot < s_min->epot) )
2553 steps_accepted++;
2555 /* Test whether the convergence criterion is met... */
2556 bDone = (s_try->fmax < inputrec->em_tol);
2558 /* Copy the arrays for force, positions and energy */
2559 /* The 'Min' array always holds the coords and forces of the minimal
2560 sampled energy */
2561 swap_em_state(&s_min, &s_try);
2562 if (count > 0)
2564 ustep *= 1.2;
2567 /* Write to trn, if necessary */
2568 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2569 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2570 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
2571 top_global, inputrec, count,
2572 s_min, state_global, energyHistory);
2574 else
2576 /* If energy is not smaller make the step smaller... */
2577 ustep *= 0.5;
2579 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2581 /* Reload the old state */
2582 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2583 s_min, top, mdatoms, fr, vsite, constr,
2584 nrnb, wcycle);
2588 /* Determine new step */
2589 stepsize = ustep/s_min->fmax;
2591 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2592 #if GMX_DOUBLE
2593 if (count == nsteps || ustep < 1e-12)
2594 #else
2595 if (count == nsteps || ustep < 1e-6)
2596 #endif
2598 if (MASTER(cr))
2600 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != nullptr);
2601 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != nullptr);
2603 bAbort = TRUE;
2606 /* Send IMD energies and positions, if bIMD is TRUE. */
2607 if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box, as_rvec_array(state_global->x.data()), inputrec, 0, wcycle) && MASTER(cr))
2609 IMD_send_positions(inputrec->imd);
2612 count++;
2613 } /* End of the loop */
2615 /* IMD cleanup, if bIMD is TRUE. */
2616 IMD_finalize(inputrec->bIMD, inputrec->imd);
2618 /* Print some data... */
2619 if (MASTER(cr))
2621 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2623 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2624 top_global, inputrec, count,
2625 s_min, state_global, energyHistory);
2627 if (MASTER(cr))
2629 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2631 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2632 s_min, sqrtNumAtoms);
2633 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2634 s_min, sqrtNumAtoms);
2637 finish_em(cr, outf, inputrec, walltime_accounting, wcycle);
2639 /* To print the actual number of steps we needed somewhere */
2640 inputrec->nsteps = count;
2642 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2644 return 0;
2645 } /* That's all folks */
2647 /*! \brief Do normal modes analysis
2648 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
2649 int nfile, const t_filenm fnm[],
2650 const gmx_output_env_t *oenv, gmx_bool bVerbose,
2651 int nstglobalcomm,
2652 gmx_vsite_t *vsite, gmx_constr_t constr,
2653 int stepout,
2654 t_inputrec *inputrec,
2655 gmx_mtop_t *top_global, t_fcdata *fcd,
2656 t_state *state_global,
2657 t_mdatoms *mdatoms,
2658 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2659 gmx_edsam_t ed,
2660 t_forcerec *fr,
2661 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
2662 real cpt_period, real max_hours,
2663 int imdport,
2664 unsigned long Flags,
2665 gmx_walltime_accounting_t walltime_accounting)
2667 double do_nm(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
2668 int nfile, const t_filenm fnm[],
2669 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
2670 int gmx_unused nstglobalcomm,
2671 gmx_vsite_t *vsite, gmx_constr_t constr,
2672 int gmx_unused stepout,
2673 t_inputrec *inputrec,
2674 gmx_mtop_t *top_global, t_fcdata *fcd,
2675 t_state *state_global,
2676 energyhistory_t gmx_unused *energyHistory,
2677 t_mdatoms *mdatoms,
2678 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2679 gmx_edsam_t gmx_unused ed,
2680 t_forcerec *fr,
2681 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2682 gmx_membed_t gmx_unused *membed,
2683 real gmx_unused cpt_period, real gmx_unused max_hours,
2684 int imdport,
2685 unsigned long gmx_unused Flags,
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);
2709 if (constr != nullptr)
2711 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2714 gmx_shellfc_t *shellfc;
2716 em_state_t state_work {};
2718 /* Init em and store the local state in state_minimum */
2719 init_em(fplog, NM, cr, inputrec,
2720 state_global, top_global, &state_work, &top,
2721 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat,
2722 vsite, constr, &shellfc,
2723 nfile, fnm, &outf, nullptr, imdport, Flags, wcycle);
2725 std::vector<size_t> atom_index = get_atom_index(top_global);
2726 snew(fneg, atom_index.size());
2727 snew(dfdx, atom_index.size());
2729 #if !GMX_DOUBLE
2730 if (bIsMaster)
2732 fprintf(stderr,
2733 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2734 " which MIGHT not be accurate enough for normal mode analysis.\n"
2735 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
2736 " are fairly modest even if you recompile in double precision.\n\n");
2738 #endif
2740 /* Check if we can/should use sparse storage format.
2742 * Sparse format is only useful when the Hessian itself is sparse, which it
2743 * will be when we use a cutoff.
2744 * For small systems (n<1000) it is easier to always use full matrix format, though.
2746 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2748 GMX_LOG(mdlog.warning).appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2749 bSparse = FALSE;
2751 else if (atom_index.size() < 1000)
2753 GMX_LOG(mdlog.warning).appendTextFormatted("Small system size (N=%d), using full Hessian format.",
2754 atom_index.size());
2755 bSparse = FALSE;
2757 else
2759 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2760 bSparse = TRUE;
2763 /* Number of dimensions, based on real atoms, that is not vsites or shell */
2764 sz = DIM*atom_index.size();
2766 fprintf(stderr, "Allocating Hessian memory...\n\n");
2768 if (bSparse)
2770 sparse_matrix = gmx_sparsematrix_init(sz);
2771 sparse_matrix->compressed_symmetric = TRUE;
2773 else
2775 snew(full_matrix, sz*sz);
2778 init_nrnb(nrnb);
2780 where();
2782 /* Write start time and temperature */
2783 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2785 /* fudge nr of steps to nr of atoms */
2786 inputrec->nsteps = atom_index.size()*2;
2788 if (bIsMaster)
2790 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2791 *(top_global->name), (int)inputrec->nsteps);
2794 nnodes = cr->nnodes;
2796 /* Make evaluate_energy do a single node force calculation */
2797 cr->nnodes = 1;
2798 evaluate_energy(fplog, cr,
2799 top_global, &state_work, top,
2800 inputrec, nrnb, wcycle, gstat,
2801 vsite, constr, fcd, graph, mdatoms, fr,
2802 mu_tot, enerd, vir, pres, -1, TRUE);
2803 cr->nnodes = nnodes;
2805 /* if forces are not small, warn user */
2806 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
2808 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
2809 if (state_work.fmax > 1.0e-3)
2811 GMX_LOG(mdlog.warning).appendText(
2812 "The force is probably not small enough to "
2813 "ensure that you are at a minimum.\n"
2814 "Be aware that negative eigenvalues may occur\n"
2815 "when the resulting matrix is diagonalized.");
2818 /***********************************************************
2820 * Loop over all pairs in matrix
2822 * do_force called twice. Once with positive and
2823 * once with negative displacement
2825 ************************************************************/
2827 /* Steps are divided one by one over the nodes */
2828 bool bNS = true;
2829 for (unsigned int aid = cr->nodeid; aid < atom_index.size(); aid += nnodes)
2831 size_t atom = atom_index[aid];
2832 for (size_t d = 0; d < DIM; d++)
2834 gmx_bool bBornRadii = FALSE;
2835 gmx_int64_t step = 0;
2836 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2837 double t = 0;
2839 x_min = state_work.s.x[atom][d];
2841 for (unsigned int dx = 0; (dx < 2); dx++)
2843 if (dx == 0)
2845 state_work.s.x[atom][d] = x_min - der_range;
2847 else
2849 state_work.s.x[atom][d] = x_min + der_range;
2852 /* Make evaluate_energy do a single node force calculation */
2853 cr->nnodes = 1;
2854 if (shellfc)
2856 /* Now is the time to relax the shells */
2857 (void) relax_shell_flexcon(fplog, cr, bVerbose, step,
2858 inputrec, bNS, force_flags,
2859 top,
2860 constr, enerd, fcd,
2861 &state_work.s, &state_work.f, vir, mdatoms,
2862 nrnb, wcycle, graph, &top_global->groups,
2863 shellfc, fr, bBornRadii, t, mu_tot,
2864 vsite);
2865 bNS = false;
2866 step++;
2868 else
2870 evaluate_energy(fplog, cr,
2871 top_global, &state_work, top,
2872 inputrec, nrnb, wcycle, gstat,
2873 vsite, constr, fcd, graph, mdatoms, fr,
2874 mu_tot, enerd, vir, pres, atom*2+dx, FALSE);
2877 cr->nnodes = nnodes;
2879 if (dx == 0)
2881 for (size_t i = 0; i < atom_index.size(); i++)
2883 copy_rvec(state_work.f[atom_index[i]], fneg[i]);
2888 /* x is restored to original */
2889 state_work.s.x[atom][d] = x_min;
2891 for (size_t j = 0; j < atom_index.size(); j++)
2893 for (size_t k = 0; (k < DIM); k++)
2895 dfdx[j][k] =
2896 -(state_work.f[atom_index[j]][k] - fneg[j][k])/(2*der_range);
2900 if (!bIsMaster)
2902 #if GMX_MPI
2903 #define mpi_type GMX_MPI_REAL
2904 MPI_Send(dfdx[0], atom_index.size()*DIM, mpi_type, MASTER(cr),
2905 cr->nodeid, cr->mpi_comm_mygroup);
2906 #endif
2908 else
2910 for (node = 0; (node < nnodes && atom+node < atom_index.size()); node++)
2912 if (node > 0)
2914 #if GMX_MPI
2915 MPI_Status stat;
2916 MPI_Recv(dfdx[0], atom_index.size()*DIM, mpi_type, node, node,
2917 cr->mpi_comm_mygroup, &stat);
2918 #undef mpi_type
2919 #endif
2922 row = (atom + node)*DIM + d;
2924 for (size_t j = 0; j < atom_index.size(); j++)
2926 for (size_t k = 0; k < DIM; k++)
2928 col = j*DIM + k;
2930 if (bSparse)
2932 if (col >= row && dfdx[j][k] != 0.0)
2934 gmx_sparsematrix_increment_value(sparse_matrix,
2935 row, col, dfdx[j][k]);
2938 else
2940 full_matrix[row*sz+col] = dfdx[j][k];
2947 if (bVerbose && fplog)
2949 fflush(fplog);
2952 /* write progress */
2953 if (bIsMaster && bVerbose)
2955 fprintf(stderr, "\rFinished step %d out of %d",
2956 static_cast<int>(std::min(atom+nnodes, atom_index.size())),
2957 static_cast<int>(atom_index.size()));
2958 fflush(stderr);
2962 if (bIsMaster)
2964 fprintf(stderr, "\n\nWriting Hessian...\n");
2965 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2968 finish_em(cr, outf, inputrec, walltime_accounting, wcycle);
2970 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size()*2);
2972 return 0;
2975 } // namespace gmx