Refactor gmx_group_t to SimulationAtomGroups
[alexxy/gromacs.git] / src / gromacs / mdrun / minimize.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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,2018,2019, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 /*! \internal \file
38  *
39  * \brief This file defines integrators for energy minimization
40  *
41  * \author Berk Hess <hess@kth.se>
42  * \author Erik Lindahl <erik@kth.se>
43  * \ingroup module_mdrun
44  */
45 #include "gmxpre.h"
46
47 #include "config.h"
48
49 #include <cmath>
50 #include <cstring>
51 #include <ctime>
52
53 #include <algorithm>
54 #include <vector>
55
56 #include "gromacs/commandline/filenm.h"
57 #include "gromacs/domdec/collect.h"
58 #include "gromacs/domdec/dlbtiming.h"
59 #include "gromacs/domdec/domdec.h"
60 #include "gromacs/domdec/domdec_struct.h"
61 #include "gromacs/domdec/partition.h"
62 #include "gromacs/ewald/pme.h"
63 #include "gromacs/fileio/confio.h"
64 #include "gromacs/fileio/mtxio.h"
65 #include "gromacs/gmxlib/network.h"
66 #include "gromacs/gmxlib/nrnb.h"
67 #include "gromacs/imd/imd.h"
68 #include "gromacs/linearalgebra/sparsematrix.h"
69 #include "gromacs/listed_forces/manage_threading.h"
70 #include "gromacs/math/functions.h"
71 #include "gromacs/math/vec.h"
72 #include "gromacs/mdlib/constr.h"
73 #include "gromacs/mdlib/dispersioncorrection.h"
74 #include "gromacs/mdlib/ebin.h"
75 #include "gromacs/mdlib/energyoutput.h"
76 #include "gromacs/mdlib/force.h"
77 #include "gromacs/mdlib/forcerec.h"
78 #include "gromacs/mdlib/gmx_omp_nthreads.h"
79 #include "gromacs/mdlib/md_support.h"
80 #include "gromacs/mdlib/mdatoms.h"
81 #include "gromacs/mdlib/mdsetup.h"
82 #include "gromacs/mdlib/ns.h"
83 #include "gromacs/mdlib/shellfc.h"
84 #include "gromacs/mdlib/stat.h"
85 #include "gromacs/mdlib/tgroup.h"
86 #include "gromacs/mdlib/trajectory_writing.h"
87 #include "gromacs/mdlib/update.h"
88 #include "gromacs/mdlib/vsite.h"
89 #include "gromacs/mdrunutility/printtime.h"
90 #include "gromacs/mdtypes/commrec.h"
91 #include "gromacs/mdtypes/inputrec.h"
92 #include "gromacs/mdtypes/md_enums.h"
93 #include "gromacs/mdtypes/mdrunoptions.h"
94 #include "gromacs/mdtypes/state.h"
95 #include "gromacs/pbcutil/mshift.h"
96 #include "gromacs/pbcutil/pbc.h"
97 #include "gromacs/timing/wallcycle.h"
98 #include "gromacs/timing/walltime_accounting.h"
99 #include "gromacs/topology/mtop_util.h"
100 #include "gromacs/topology/topology.h"
101 #include "gromacs/utility/cstringutil.h"
102 #include "gromacs/utility/exceptions.h"
103 #include "gromacs/utility/fatalerror.h"
104 #include "gromacs/utility/logger.h"
105 #include "gromacs/utility/smalloc.h"
106
107 #include "integrator.h"
108
109 //! Utility structure for manipulating states during EM
110 typedef struct {
111     //! Copy of the global state
112     t_state                 s;
113     //! Force array
114     PaddedVector<gmx::RVec> f;
115     //! Potential energy
116     real                    epot;
117     //! Norm of the force
118     real                    fnorm;
119     //! Maximum force
120     real                    fmax;
121     //! Direction
122     int                     a_fmax;
123 } em_state_t;
124
125 //! Print the EM starting conditions
126 static void print_em_start(FILE                     *fplog,
127                            const t_commrec          *cr,
128                            gmx_walltime_accounting_t walltime_accounting,
129                            gmx_wallcycle_t           wcycle,
130                            const char               *name)
131 {
132     walltime_accounting_start_time(walltime_accounting);
133     wallcycle_start(wcycle, ewcRUN);
134     print_start(fplog, cr, walltime_accounting, name);
135 }
136
137 //! Stop counting time for EM
138 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
139                         gmx_wallcycle_t           wcycle)
140 {
141     wallcycle_stop(wcycle, ewcRUN);
142
143     walltime_accounting_end_time(walltime_accounting);
144 }
145
146 //! Printing a log file and console header
147 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
148 {
149     fprintf(out, "\n");
150     fprintf(out, "%s:\n", minimizer);
151     fprintf(out, "   Tolerance (Fmax)   = %12.5e\n", ftol);
152     fprintf(out, "   Number of steps    = %12d\n", nsteps);
153 }
154
155 //! Print warning message
156 static void warn_step(FILE     *fp,
157                       real      ftol,
158                       real      fmax,
159                       gmx_bool  bLastStep,
160                       gmx_bool  bConstrain)
161 {
162     constexpr bool realIsDouble = GMX_DOUBLE;
163     char           buffer[2048];
164
165     if (!std::isfinite(fmax))
166     {
167         sprintf(buffer,
168                 "\nEnergy minimization has stopped because the force "
169                 "on at least one atom is not finite. This usually means "
170                 "atoms are overlapping. Modify the input coordinates to "
171                 "remove atom overlap or use soft-core potentials with "
172                 "the free energy code to avoid infinite forces.\n%s",
173                 !realIsDouble ?
174                 "You could also be lucky that switching to double precision "
175                 "is sufficient to obtain finite forces.\n" :
176                 "");
177     }
178     else if (bLastStep)
179     {
180         sprintf(buffer,
181                 "\nEnergy minimization reached the maximum number "
182                 "of steps before the forces reached the requested "
183                 "precision Fmax < %g.\n", ftol);
184     }
185     else
186     {
187         sprintf(buffer,
188                 "\nEnergy minimization has stopped, but the forces have "
189                 "not converged to the requested precision Fmax < %g (which "
190                 "may not be possible for your system). It stopped "
191                 "because the algorithm tried to make a new step whose size "
192                 "was too small, or there was no change in the energy since "
193                 "last step. Either way, we regard the minimization as "
194                 "converged to within the available machine precision, "
195                 "given your starting configuration and EM parameters.\n%s%s",
196                 ftol,
197                 !realIsDouble ?
198                 "\nDouble precision normally gives you higher accuracy, but "
199                 "this is often not needed for preparing to run molecular "
200                 "dynamics.\n" :
201                 "",
202                 bConstrain ?
203                 "You might need to increase your constraint accuracy, or turn\n"
204                 "off constraints altogether (set constraints = none in mdp file)\n" :
205                 "");
206     }
207
208     fputs(wrap_lines(buffer, 78, 0, FALSE), stderr);
209     fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
210 }
211
212 //! Print message about convergence of the EM
213 static void print_converged(FILE *fp, const char *alg, real ftol,
214                             int64_t count, gmx_bool bDone, int64_t nsteps,
215                             const em_state_t *ems, double sqrtNumAtoms)
216 {
217     char buf[STEPSTRSIZE];
218
219     if (bDone)
220     {
221         fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
222                 alg, ftol, gmx_step_str(count, buf));
223     }
224     else if (count < nsteps)
225     {
226         fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
227                 "but did not reach the requested Fmax < %g.\n",
228                 alg, gmx_step_str(count, buf), ftol);
229     }
230     else
231     {
232         fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
233                 alg, ftol, gmx_step_str(count, buf));
234     }
235
236 #if GMX_DOUBLE
237     fprintf(fp, "Potential Energy  = %21.14e\n", ems->epot);
238     fprintf(fp, "Maximum force     = %21.14e on atom %d\n", ems->fmax, ems->a_fmax + 1);
239     fprintf(fp, "Norm of force     = %21.14e\n", ems->fnorm/sqrtNumAtoms);
240 #else
241     fprintf(fp, "Potential Energy  = %14.7e\n", ems->epot);
242     fprintf(fp, "Maximum force     = %14.7e on atom %d\n", ems->fmax, ems->a_fmax + 1);
243     fprintf(fp, "Norm of force     = %14.7e\n", ems->fnorm/sqrtNumAtoms);
244 #endif
245 }
246
247 //! Compute the norm and max of the force array in parallel
248 static void get_f_norm_max(const t_commrec *cr,
249                            t_grpopts *opts, t_mdatoms *mdatoms, const rvec *f,
250                            real *fnorm, real *fmax, int *a_fmax)
251 {
252     double fnorm2, *sum;
253     real   fmax2, fam;
254     int    la_max, a_max, start, end, i, m, gf;
255
256     /* This routine finds the largest force and returns it.
257      * On parallel machines the global max is taken.
258      */
259     fnorm2 = 0;
260     fmax2  = 0;
261     la_max = -1;
262     start  = 0;
263     end    = mdatoms->homenr;
264     if (mdatoms->cFREEZE)
265     {
266         for (i = start; i < end; i++)
267         {
268             gf  = mdatoms->cFREEZE[i];
269             fam = 0;
270             for (m = 0; m < DIM; m++)
271             {
272                 if (!opts->nFreeze[gf][m])
273                 {
274                     fam += gmx::square(f[i][m]);
275                 }
276             }
277             fnorm2 += fam;
278             if (fam > fmax2)
279             {
280                 fmax2  = fam;
281                 la_max = i;
282             }
283         }
284     }
285     else
286     {
287         for (i = start; i < end; i++)
288         {
289             fam     = norm2(f[i]);
290             fnorm2 += fam;
291             if (fam > fmax2)
292             {
293                 fmax2  = fam;
294                 la_max = i;
295             }
296         }
297     }
298
299     if (la_max >= 0 && DOMAINDECOMP(cr))
300     {
301         a_max = cr->dd->globalAtomIndices[la_max];
302     }
303     else
304     {
305         a_max = la_max;
306     }
307     if (PAR(cr))
308     {
309         snew(sum, 2*cr->nnodes+1);
310         sum[2*cr->nodeid]   = fmax2;
311         sum[2*cr->nodeid+1] = a_max;
312         sum[2*cr->nnodes]   = fnorm2;
313         gmx_sumd(2*cr->nnodes+1, sum, cr);
314         fnorm2 = sum[2*cr->nnodes];
315         /* Determine the global maximum */
316         for (i = 0; i < cr->nnodes; i++)
317         {
318             if (sum[2*i] > fmax2)
319             {
320                 fmax2 = sum[2*i];
321                 a_max = gmx::roundToInt(sum[2*i+1]);
322             }
323         }
324         sfree(sum);
325     }
326
327     if (fnorm)
328     {
329         *fnorm = sqrt(fnorm2);
330     }
331     if (fmax)
332     {
333         *fmax  = sqrt(fmax2);
334     }
335     if (a_fmax)
336     {
337         *a_fmax = a_max;
338     }
339 }
340
341 //! Compute the norm of the force
342 static void get_state_f_norm_max(const t_commrec *cr,
343                                  t_grpopts *opts, t_mdatoms *mdatoms,
344                                  em_state_t *ems)
345 {
346     get_f_norm_max(cr, opts, mdatoms, ems->f.rvec_array(),
347                    &ems->fnorm, &ems->fmax, &ems->a_fmax);
348 }
349
350 //! Initialize the energy minimization
351 static void init_em(FILE *fplog,
352                     const gmx::MDLogger &mdlog,
353                     const char *title,
354                     const t_commrec *cr,
355                     const gmx_multisim_t *ms,
356                     t_inputrec *ir,
357                     const gmx::MdrunOptions &mdrunOptions,
358                     t_state *state_global, gmx_mtop_t *top_global,
359                     em_state_t *ems, gmx_localtop_t *top,
360                     t_nrnb *nrnb,
361                     t_forcerec *fr,
362                     t_graph **graph, gmx::MDAtoms *mdAtoms, gmx_global_stat_t *gstat,
363                     gmx_vsite_t *vsite, gmx::Constraints *constr, gmx_shellfc_t **shellfc,
364                     int nfile, const t_filenm fnm[])
365 {
366     real dvdl_constr;
367
368     if (fplog)
369     {
370         fprintf(fplog, "Initiating %s\n", title);
371     }
372
373     if (MASTER(cr))
374     {
375         state_global->ngtc = 0;
376     }
377     initialize_lambdas(fplog, *ir, MASTER(cr), &(state_global->fep_state), state_global->lambda, nullptr);
378
379     init_nrnb(nrnb);
380
381     /* Interactive molecular dynamics */
382     init_IMD(ir, cr, ms, top_global, fplog, 1,
383              MASTER(cr) ? state_global->x.rvec_array() : nullptr,
384              nfile, fnm, nullptr, mdrunOptions);
385
386     if (ir->eI == eiNM)
387     {
388         GMX_ASSERT(shellfc != nullptr, "With NM we always support shells");
389
390         *shellfc = init_shell_flexcon(stdout,
391                                       top_global,
392                                       constr ? constr->numFlexibleConstraints() : 0,
393                                       ir->nstcalcenergy,
394                                       DOMAINDECOMP(cr));
395     }
396     else
397     {
398         GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI), "This else currently only handles energy minimizers, consider if your algorithm needs shell/flexible-constraint support");
399
400         /* With energy minimization, shells and flexible constraints are
401          * automatically minimized when treated like normal DOFS.
402          */
403         if (shellfc != nullptr)
404         {
405             *shellfc = nullptr;
406         }
407     }
408
409     auto mdatoms = mdAtoms->mdatoms();
410     if (DOMAINDECOMP(cr))
411     {
412         top->useInDomainDecomp_ = true;
413         dd_init_local_top(*top_global, top);
414
415         dd_init_local_state(cr->dd, state_global, &ems->s);
416
417         /* Distribute the charge groups over the nodes from the master node */
418         dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1,
419                             state_global, *top_global, ir,
420                             &ems->s, &ems->f, mdAtoms, top,
421                             fr, vsite, constr,
422                             nrnb, nullptr, FALSE);
423         dd_store_state(cr->dd, &ems->s);
424
425         *graph = nullptr;
426     }
427     else
428     {
429         state_change_natoms(state_global, state_global->natoms);
430         /* Just copy the state */
431         ems->s = *state_global;
432         state_change_natoms(&ems->s, ems->s.natoms);
433         ems->f.resizeWithPadding(ems->s.natoms);
434
435         mdAlgorithmsSetupAtomData(cr, ir, *top_global, top, fr,
436                                   graph, mdAtoms,
437                                   constr, vsite, shellfc ? *shellfc : nullptr);
438
439         if (vsite)
440         {
441             set_vsite_top(vsite, top, mdatoms);
442         }
443     }
444
445     update_mdatoms(mdAtoms->mdatoms(), ems->s.lambda[efptMASS]);
446
447     if (constr)
448     {
449         // TODO how should this cross-module support dependency be managed?
450         if (ir->eConstrAlg == econtSHAKE &&
451             gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
452         {
453             gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
454                       econstr_names[econtSHAKE], econstr_names[econtLINCS]);
455         }
456
457         if (!ir->bContinuation)
458         {
459             /* Constrain the starting coordinates */
460             dvdl_constr = 0;
461             constr->apply(TRUE, TRUE,
462                           -1, 0, 1.0,
463                           ems->s.x.rvec_array(),
464                           ems->s.x.rvec_array(),
465                           nullptr,
466                           ems->s.box,
467                           ems->s.lambda[efptFEP], &dvdl_constr,
468                           nullptr, nullptr, gmx::ConstraintVariable::Positions);
469         }
470     }
471
472     if (PAR(cr))
473     {
474         *gstat = global_stat_init(ir);
475     }
476     else
477     {
478         *gstat = nullptr;
479     }
480
481     calc_shifts(ems->s.box, fr->shift_vec);
482 }
483
484 //! Finalize the minimization
485 static void finish_em(const t_commrec *cr, gmx_mdoutf_t outf,
486                       gmx_walltime_accounting_t walltime_accounting,
487                       gmx_wallcycle_t wcycle)
488 {
489     if (!thisRankHasDuty(cr, DUTY_PME))
490     {
491         /* Tell the PME only node to finish */
492         gmx_pme_send_finish(cr);
493     }
494
495     done_mdoutf(outf);
496
497     em_time_end(walltime_accounting, wcycle);
498 }
499
500 //! Swap two different EM states during minimization
501 static void swap_em_state(em_state_t **ems1, em_state_t **ems2)
502 {
503     em_state_t *tmp;
504
505     tmp   = *ems1;
506     *ems1 = *ems2;
507     *ems2 = tmp;
508 }
509
510 //! Save the EM trajectory
511 static void write_em_traj(FILE *fplog, const t_commrec *cr,
512                           gmx_mdoutf_t outf,
513                           gmx_bool bX, gmx_bool bF, const char *confout,
514                           gmx_mtop_t *top_global,
515                           t_inputrec *ir, int64_t step,
516                           em_state_t *state,
517                           t_state *state_global,
518                           ObservablesHistory *observablesHistory)
519 {
520     int mdof_flags = 0;
521
522     if (bX)
523     {
524         mdof_flags |= MDOF_X;
525     }
526     if (bF)
527     {
528         mdof_flags |= MDOF_F;
529     }
530
531     /* If we want IMD output, set appropriate MDOF flag */
532     if (ir->bIMD)
533     {
534         mdof_flags |= MDOF_IMD;
535     }
536
537     mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
538                                      top_global, step, static_cast<double>(step),
539                                      &state->s, state_global, observablesHistory,
540                                      state->f);
541
542     if (confout != nullptr)
543     {
544         if (DOMAINDECOMP(cr))
545         {
546             /* If bX=true, x was collected to state_global in the call above */
547             if (!bX)
548             {
549                 auto globalXRef = MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>();
550                 dd_collect_vec(cr->dd, &state->s, state->s.x, globalXRef);
551             }
552         }
553         else
554         {
555             /* Copy the local state pointer */
556             state_global = &state->s;
557         }
558
559         if (MASTER(cr))
560         {
561             if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
562             {
563                 /* Make molecules whole only for confout writing */
564                 do_pbc_mtop(ir->ePBC, state->s.box, top_global,
565                             state_global->x.rvec_array());
566             }
567
568             write_sto_conf_mtop(confout,
569                                 *top_global->name, top_global,
570                                 state_global->x.rvec_array(), nullptr, ir->ePBC, state->s.box);
571         }
572     }
573 }
574
575 //! \brief Do one minimization step
576 //
577 // \returns true when the step succeeded, false when a constraint error occurred
578 static bool do_em_step(const t_commrec *cr,
579                        t_inputrec *ir, t_mdatoms *md,
580                        em_state_t *ems1, real a, const PaddedVector<gmx::RVec> *force,
581                        em_state_t *ems2,
582                        gmx::Constraints *constr,
583                        int64_t count)
584
585 {
586     t_state *s1, *s2;
587     int      start, end;
588     real     dvdl_constr;
589     int      nthreads gmx_unused;
590
591     bool     validStep = true;
592
593     s1 = &ems1->s;
594     s2 = &ems2->s;
595
596     if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
597     {
598         gmx_incons("state mismatch in do_em_step");
599     }
600
601     s2->flags = s1->flags;
602
603     if (s2->natoms != s1->natoms)
604     {
605         state_change_natoms(s2, s1->natoms);
606         ems2->f.resizeWithPadding(s2->natoms);
607     }
608     if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
609     {
610         s2->cg_gl.resize(s1->cg_gl.size());
611     }
612
613     copy_mat(s1->box, s2->box);
614     /* Copy free energy state */
615     s2->lambda = s1->lambda;
616     copy_mat(s1->box, s2->box);
617
618     start = 0;
619     end   = md->homenr;
620
621     nthreads = gmx_omp_nthreads_get(emntUpdate);
622 #pragma omp parallel num_threads(nthreads)
623     {
624         const rvec *x1 = s1->x.rvec_array();
625         rvec       *x2 = s2->x.rvec_array();
626         const rvec *f  = force->rvec_array();
627
628         int         gf = 0;
629 #pragma omp for schedule(static) nowait
630         for (int i = start; i < end; i++)
631         {
632             try
633             {
634                 if (md->cFREEZE)
635                 {
636                     gf = md->cFREEZE[i];
637                 }
638                 for (int m = 0; m < DIM; m++)
639                 {
640                     if (ir->opts.nFreeze[gf][m])
641                     {
642                         x2[i][m] = x1[i][m];
643                     }
644                     else
645                     {
646                         x2[i][m] = x1[i][m] + a*f[i][m];
647                     }
648                 }
649             }
650             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
651         }
652
653         if (s2->flags & (1<<estCGP))
654         {
655             /* Copy the CG p vector */
656             const rvec *p1 = s1->cg_p.rvec_array();
657             rvec       *p2 = s2->cg_p.rvec_array();
658 #pragma omp for schedule(static) nowait
659             for (int i = start; i < end; i++)
660             {
661                 // Trivial OpenMP block that does not throw
662                 copy_rvec(p1[i], p2[i]);
663             }
664         }
665
666         if (DOMAINDECOMP(cr))
667         {
668             s2->ddp_count = s1->ddp_count;
669
670             /* OpenMP does not supported unsigned loop variables */
671 #pragma omp for schedule(static) nowait
672             for (int i = 0; i < gmx::ssize(s2->cg_gl); i++)
673             {
674                 s2->cg_gl[i] = s1->cg_gl[i];
675             }
676             s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
677         }
678     }
679
680     if (constr)
681     {
682         dvdl_constr = 0;
683         validStep   =
684             constr->apply(TRUE, TRUE,
685                           count, 0, 1.0,
686                           s1->x.rvec_array(), s2->x.rvec_array(),
687                           nullptr, s2->box,
688                           s2->lambda[efptBONDED], &dvdl_constr,
689                           nullptr, nullptr, gmx::ConstraintVariable::Positions);
690
691         if (cr->nnodes > 1)
692         {
693             /* This global reduction will affect performance at high
694              * parallelization, but we can not really avoid it.
695              * But usually EM is not run at high parallelization.
696              */
697             int reductionBuffer = static_cast<int>(!validStep);
698             gmx_sumi(1, &reductionBuffer, cr);
699             validStep           = (reductionBuffer == 0);
700         }
701
702         // We should move this check to the different minimizers
703         if (!validStep && ir->eI != eiSteep)
704         {
705             gmx_fatal(FARGS, "The coordinates could not be constrained. Minimizer '%s' can not handle constraint failures, use minimizer '%s' before using '%s'.",
706                       EI(ir->eI), EI(eiSteep), EI(ir->eI));
707         }
708     }
709
710     return validStep;
711 }
712
713 //! Prepare EM for using domain decomposition parallellization
714 static void em_dd_partition_system(FILE *fplog,
715                                    const gmx::MDLogger &mdlog,
716                                    int step, const t_commrec *cr,
717                                    gmx_mtop_t *top_global, t_inputrec *ir,
718                                    em_state_t *ems, gmx_localtop_t *top,
719                                    gmx::MDAtoms *mdAtoms, t_forcerec *fr,
720                                    gmx_vsite_t *vsite, gmx::Constraints *constr,
721                                    t_nrnb *nrnb, gmx_wallcycle_t wcycle)
722 {
723     /* Repartition the domain decomposition */
724     dd_partition_system(fplog, mdlog, step, cr, FALSE, 1,
725                         nullptr, *top_global, ir,
726                         &ems->s, &ems->f,
727                         mdAtoms, top, fr, vsite, constr,
728                         nrnb, wcycle, FALSE);
729     dd_store_state(cr->dd, &ems->s);
730 }
731
732 namespace
733 {
734
735 /*! \brief Class to handle the work of setting and doing an energy evaluation.
736  *
737  * This class is a mere aggregate of parameters to pass to evaluate an
738  * energy, so that future changes to names and types of them consume
739  * less time when refactoring other code.
740  *
741  * Aggregate initialization is used, for which the chief risk is that
742  * if a member is added at the end and not all initializer lists are
743  * updated, then the member will be value initialized, which will
744  * typically mean initialization to zero.
745  *
746  * We only want to construct one of these with an initializer list, so
747  * we explicitly delete the default constructor. */
748 class EnergyEvaluator
749 {
750     public:
751         //! We only intend to construct such objects with an initializer list.
752         EnergyEvaluator() = delete;
753         /*! \brief Evaluates an energy on the state in \c ems.
754          *
755          * \todo In practice, the same objects mu_tot, vir, and pres
756          * are always passed to this function, so we would rather have
757          * them as data members. However, their C-array types are
758          * unsuited for aggregate initialization. When the types
759          * improve, the call signature of this method can be reduced.
760          */
761         void run(em_state_t *ems, rvec mu_tot,
762                  tensor vir, tensor pres,
763                  int64_t count, gmx_bool bFirst);
764         //! Handles logging (deprecated).
765         FILE                 *fplog;
766         //! Handles logging.
767         const gmx::MDLogger  &mdlog;
768         //! Handles communication.
769         const t_commrec      *cr;
770         //! Coordinates multi-simulations.
771         const gmx_multisim_t *ms;
772         //! Holds the simulation topology.
773         gmx_mtop_t           *top_global;
774         //! Holds the domain topology.
775         gmx_localtop_t       *top;
776         //! User input options.
777         t_inputrec           *inputrec;
778         //! Manages flop accounting.
779         t_nrnb               *nrnb;
780         //! Manages wall cycle accounting.
781         gmx_wallcycle_t       wcycle;
782         //! Coordinates global reduction.
783         gmx_global_stat_t     gstat;
784         //! Handles virtual sites.
785         gmx_vsite_t          *vsite;
786         //! Handles constraints.
787         gmx::Constraints     *constr;
788         //! Handles strange things.
789         t_fcdata             *fcd;
790         //! Molecular graph for SHAKE.
791         t_graph              *graph;
792         //! Per-atom data for this domain.
793         gmx::MDAtoms         *mdAtoms;
794         //! Handles how to calculate the forces.
795         t_forcerec           *fr;
796         //! Schedule of force-calculation work each step for this task.
797         gmx::PpForceWorkload *ppForceWorkload;
798         //! Stores the computed energies.
799         gmx_enerdata_t       *enerd;
800 };
801
802 void
803 EnergyEvaluator::run(em_state_t *ems, rvec mu_tot,
804                      tensor vir, tensor pres,
805                      int64_t count, gmx_bool bFirst)
806 {
807     real     t;
808     gmx_bool bNS;
809     tensor   force_vir, shake_vir, ekin;
810     real     dvdl_constr, prescorr, enercorr, dvdlcorr;
811     real     terminate = 0;
812
813     /* Set the time to the initial time, the time does not change during EM */
814     t = inputrec->init_t;
815
816     if (bFirst ||
817         (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
818     {
819         /* This is the first state or an old state used before the last ns */
820         bNS = TRUE;
821     }
822     else
823     {
824         bNS = FALSE;
825         if (inputrec->nstlist > 0)
826         {
827             bNS = TRUE;
828         }
829     }
830
831     if (vsite)
832     {
833         construct_vsites(vsite, ems->s.x.rvec_array(), 1, nullptr,
834                          top->idef.iparams, top->idef.il,
835                          fr->ePBC, fr->bMolPBC, cr, ems->s.box);
836     }
837
838     if (DOMAINDECOMP(cr) && bNS)
839     {
840         /* Repartition the domain decomposition */
841         em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec,
842                                ems, top, mdAtoms, fr, vsite, constr,
843                                nrnb, wcycle);
844     }
845
846     /* Calc force & energy on new trial position  */
847     /* do_force always puts the charge groups in the box and shifts again
848      * We do not unshift, so molecules are always whole in congrad.c
849      */
850     do_force(fplog, cr, ms, inputrec, nullptr, nullptr,
851              count, nrnb, wcycle, top, &top_global->groups,
852              ems->s.box, ems->s.x.arrayRefWithPadding(), &ems->s.hist,
853              ems->f.arrayRefWithPadding(), force_vir, mdAtoms->mdatoms(), enerd, fcd,
854              ems->s.lambda, graph, fr, ppForceWorkload, vsite, mu_tot, t, nullptr,
855              GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
856              GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
857              (bNS ? GMX_FORCE_NS : 0),
858              DDBalanceRegionHandler(cr));
859
860     /* Clear the unused shake virial and pressure */
861     clear_mat(shake_vir);
862     clear_mat(pres);
863
864     /* Communicate stuff when parallel */
865     if (PAR(cr) && inputrec->eI != eiNM)
866     {
867         wallcycle_start(wcycle, ewcMoveE);
868
869         global_stat(gstat, cr, enerd, force_vir, shake_vir, mu_tot,
870                     inputrec, nullptr, nullptr, nullptr, 1, &terminate,
871                     nullptr, FALSE,
872                     CGLO_ENERGY |
873                     CGLO_PRESSURE |
874                     CGLO_CONSTRAINT);
875
876         wallcycle_stop(wcycle, ewcMoveE);
877     }
878
879     /* Calculate long range corrections to pressure and energy */
880     calc_dispcorr(inputrec, fr, ems->s.box, ems->s.lambda[efptVDW],
881                   pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
882     enerd->term[F_DISPCORR] = enercorr;
883     enerd->term[F_EPOT]    += enercorr;
884     enerd->term[F_PRES]    += prescorr;
885     enerd->term[F_DVDL]    += dvdlcorr;
886
887     ems->epot = enerd->term[F_EPOT];
888
889     if (constr)
890     {
891         /* Project out the constraint components of the force */
892         dvdl_constr = 0;
893         rvec *f_rvec = ems->f.rvec_array();
894         constr->apply(FALSE, FALSE,
895                       count, 0, 1.0,
896                       ems->s.x.rvec_array(), f_rvec, f_rvec,
897                       ems->s.box,
898                       ems->s.lambda[efptBONDED], &dvdl_constr,
899                       nullptr, &shake_vir, gmx::ConstraintVariable::ForceDispl);
900         enerd->term[F_DVDL_CONSTR] += dvdl_constr;
901         m_add(force_vir, shake_vir, vir);
902     }
903     else
904     {
905         copy_mat(force_vir, vir);
906     }
907
908     clear_mat(ekin);
909     enerd->term[F_PRES] =
910         calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
911
912     sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
913
914     if (EI_ENERGY_MINIMIZATION(inputrec->eI))
915     {
916         get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
917     }
918 }
919
920 } // namespace
921
922 //! Parallel utility summing energies and forces
923 static double reorder_partsum(const t_commrec *cr, t_grpopts *opts,
924                               gmx_mtop_t *top_global,
925                               em_state_t *s_min, em_state_t *s_b)
926 {
927     t_block       *cgs_gl;
928     int            ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
929     double         partsum;
930
931     if (debug)
932     {
933         fprintf(debug, "Doing reorder_partsum\n");
934     }
935
936     const rvec *fm = s_min->f.rvec_array();
937     const rvec *fb = s_b->f.rvec_array();
938
939     cgs_gl = dd_charge_groups_global(cr->dd);
940     index  = cgs_gl->index;
941
942     /* Collect fm in a global vector fmg.
943      * This conflicts with the spirit of domain decomposition,
944      * but to fully optimize this a much more complicated algorithm is required.
945      */
946     rvec *fmg;
947     snew(fmg, top_global->natoms);
948
949     ncg   = s_min->s.cg_gl.size();
950     cg_gl = s_min->s.cg_gl.data();
951     i     = 0;
952     for (c = 0; c < ncg; c++)
953     {
954         cg = cg_gl[c];
955         a0 = index[cg];
956         a1 = index[cg+1];
957         for (a = a0; a < a1; a++)
958         {
959             copy_rvec(fm[i], fmg[a]);
960             i++;
961         }
962     }
963     gmx_sum(top_global->natoms*3, fmg[0], cr);
964
965     /* Now we will determine the part of the sum for the cgs in state s_b */
966     ncg         = s_b->s.cg_gl.size();
967     cg_gl       = s_b->s.cg_gl.data();
968     partsum     = 0;
969     i           = 0;
970     gf          = 0;
971     gmx::ArrayRef<unsigned char> grpnrFREEZE = top_global->groups.groupNumbers[SimulationAtomGroupType::Freeze];
972     for (c = 0; c < ncg; c++)
973     {
974         cg = cg_gl[c];
975         a0 = index[cg];
976         a1 = index[cg+1];
977         for (a = a0; a < a1; a++)
978         {
979             if (!grpnrFREEZE.empty())
980             {
981                 gf = grpnrFREEZE[i];
982             }
983             for (m = 0; m < DIM; m++)
984             {
985                 if (!opts->nFreeze[gf][m])
986                 {
987                     partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
988                 }
989             }
990             i++;
991         }
992     }
993
994     sfree(fmg);
995
996     return partsum;
997 }
998
999 //! Print some stuff, like beta, whatever that means.
1000 static real pr_beta(const t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
1001                     gmx_mtop_t *top_global,
1002                     em_state_t *s_min, em_state_t *s_b)
1003 {
1004     double sum;
1005
1006     /* This is just the classical Polak-Ribiere calculation of beta;
1007      * it looks a bit complicated since we take freeze groups into account,
1008      * and might have to sum it in parallel runs.
1009      */
1010
1011     if (!DOMAINDECOMP(cr) ||
1012         (s_min->s.ddp_count == cr->dd->ddp_count &&
1013          s_b->s.ddp_count   == cr->dd->ddp_count))
1014     {
1015         const rvec *fm  = s_min->f.rvec_array();
1016         const rvec *fb  = s_b->f.rvec_array();
1017         sum             = 0;
1018         int         gf  = 0;
1019         /* This part of code can be incorrect with DD,
1020          * since the atom ordering in s_b and s_min might differ.
1021          */
1022         for (int i = 0; i < mdatoms->homenr; i++)
1023         {
1024             if (mdatoms->cFREEZE)
1025             {
1026                 gf = mdatoms->cFREEZE[i];
1027             }
1028             for (int m = 0; m < DIM; m++)
1029             {
1030                 if (!opts->nFreeze[gf][m])
1031                 {
1032                     sum += (fb[i][m] - fm[i][m])*fb[i][m];
1033                 }
1034             }
1035         }
1036     }
1037     else
1038     {
1039         /* We need to reorder cgs while summing */
1040         sum = reorder_partsum(cr, opts, top_global, s_min, s_b);
1041     }
1042     if (PAR(cr))
1043     {
1044         gmx_sumd(1, &sum, cr);
1045     }
1046
1047     return sum/gmx::square(s_min->fnorm);
1048 }
1049
1050 namespace gmx
1051 {
1052
1053 void
1054 Integrator::do_cg()
1055 {
1056     const char        *CG = "Polak-Ribiere Conjugate Gradients";
1057
1058     gmx_localtop_t     top;
1059     gmx_global_stat_t  gstat;
1060     t_graph           *graph;
1061     double             tmp, minstep;
1062     real               stepsize;
1063     real               a, b, c, beta = 0.0;
1064     real               epot_repl = 0;
1065     real               pnorm;
1066     gmx_bool           converged, foundlower;
1067     rvec               mu_tot = {0};
1068     gmx_bool           do_log = FALSE, do_ene = FALSE, do_x, do_f;
1069     tensor             vir, pres;
1070     int                number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1071     int                m, step, nminstep;
1072     auto               mdatoms = mdAtoms->mdatoms();
1073
1074     GMX_LOG(mdlog.info).asParagraph().
1075         appendText("Note that activating conjugate gradient energy minimization via the "
1076                    "integrator .mdp option and the command gmx mdrun may "
1077                    "be available in a different form in a future version of GROMACS, "
1078                    "e.g. gmx minimize and an .mdp option.");
1079
1080     step = 0;
1081
1082     if (MASTER(cr))
1083     {
1084         // In CG, the state is extended with a search direction
1085         state_global->flags |= (1<<estCGP);
1086
1087         // Ensure the extra per-atom state array gets allocated
1088         state_change_natoms(state_global, state_global->natoms);
1089
1090         // Initialize the search direction to zero
1091         for (RVec &cg_p : state_global->cg_p)
1092         {
1093             cg_p = { 0, 0, 0 };
1094         }
1095     }
1096
1097     /* Create 4 states on the stack and extract pointers that we will swap */
1098     em_state_t  s0 {}, s1 {}, s2 {}, s3 {};
1099     em_state_t *s_min = &s0;
1100     em_state_t *s_a   = &s1;
1101     em_state_t *s_b   = &s2;
1102     em_state_t *s_c   = &s3;
1103
1104     /* Init em and store the local state in s_min */
1105     init_em(fplog, mdlog, CG, cr, ms, inputrec, mdrunOptions,
1106             state_global, top_global, s_min, &top,
1107             nrnb, fr, &graph, mdAtoms, &gstat,
1108             vsite, constr, nullptr,
1109             nfile, fnm);
1110     gmx_mdoutf       *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, inputrec, top_global, nullptr, wcycle);
1111     gmx::EnergyOutput energyOutput;
1112     energyOutput.prepare(mdoutf_get_fp_ene(outf), top_global, inputrec, nullptr);
1113
1114     /* Print to log file */
1115     print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1116
1117     /* Max number of steps */
1118     number_steps = inputrec->nsteps;
1119
1120     if (MASTER(cr))
1121     {
1122         sp_header(stderr, CG, inputrec->em_tol, number_steps);
1123     }
1124     if (fplog)
1125     {
1126         sp_header(fplog, CG, inputrec->em_tol, number_steps);
1127     }
1128
1129     EnergyEvaluator energyEvaluator {
1130         fplog, mdlog, cr, ms,
1131         top_global, &top,
1132         inputrec, nrnb, wcycle, gstat,
1133         vsite, constr, fcd, graph,
1134         mdAtoms, fr, ppForceWorkload, enerd
1135     };
1136     /* Call the force routine and some auxiliary (neighboursearching etc.) */
1137     /* do_force always puts the charge groups in the box and shifts again
1138      * We do not unshift, so molecules are always whole in congrad.c
1139      */
1140     energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE);
1141
1142     if (MASTER(cr))
1143     {
1144         /* Copy stuff to the energy bin for easy printing etc. */
1145         matrix nullBox = {};
1146         energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step),
1147                                          mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1148                                          nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1149
1150         print_ebin_header(fplog, step, step);
1151         energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE,
1152                                            fplog, step, step, eprNORMAL,
1153                                            fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1154     }
1155
1156     /* Estimate/guess the initial stepsize */
1157     stepsize = inputrec->em_stepsize/s_min->fnorm;
1158
1159     if (MASTER(cr))
1160     {
1161         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1162         fprintf(stderr, "   F-max             = %12.5e on atom %d\n",
1163                 s_min->fmax, s_min->a_fmax+1);
1164         fprintf(stderr, "   F-Norm            = %12.5e\n",
1165                 s_min->fnorm/sqrtNumAtoms);
1166         fprintf(stderr, "\n");
1167         /* and copy to the log file too... */
1168         fprintf(fplog, "   F-max             = %12.5e on atom %d\n",
1169                 s_min->fmax, s_min->a_fmax+1);
1170         fprintf(fplog, "   F-Norm            = %12.5e\n",
1171                 s_min->fnorm/sqrtNumAtoms);
1172         fprintf(fplog, "\n");
1173     }
1174     /* Start the loop over CG steps.
1175      * Each successful step is counted, and we continue until
1176      * we either converge or reach the max number of steps.
1177      */
1178     converged = FALSE;
1179     for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1180     {
1181
1182         /* start taking steps in a new direction
1183          * First time we enter the routine, beta=0, and the direction is
1184          * simply the negative gradient.
1185          */
1186
1187         /* Calculate the new direction in p, and the gradient in this direction, gpa */
1188         rvec       *pm  = s_min->s.cg_p.rvec_array();
1189         const rvec *sfm = s_min->f.rvec_array();
1190         double      gpa = 0;
1191         int         gf  = 0;
1192         for (int i = 0; i < mdatoms->homenr; i++)
1193         {
1194             if (mdatoms->cFREEZE)
1195             {
1196                 gf = mdatoms->cFREEZE[i];
1197             }
1198             for (m = 0; m < DIM; m++)
1199             {
1200                 if (!inputrec->opts.nFreeze[gf][m])
1201                 {
1202                     pm[i][m] = sfm[i][m] + beta*pm[i][m];
1203                     gpa     -= pm[i][m]*sfm[i][m];
1204                     /* f is negative gradient, thus the sign */
1205                 }
1206                 else
1207                 {
1208                     pm[i][m] = 0;
1209                 }
1210             }
1211         }
1212
1213         /* Sum the gradient along the line across CPUs */
1214         if (PAR(cr))
1215         {
1216             gmx_sumd(1, &gpa, cr);
1217         }
1218
1219         /* Calculate the norm of the search vector */
1220         get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1221
1222         /* Just in case stepsize reaches zero due to numerical precision... */
1223         if (stepsize <= 0)
1224         {
1225             stepsize = inputrec->em_stepsize/pnorm;
1226         }
1227
1228         /*
1229          * Double check the value of the derivative in the search direction.
1230          * If it is positive it must be due to the old information in the
1231          * CG formula, so just remove that and start over with beta=0.
1232          * This corresponds to a steepest descent step.
1233          */
1234         if (gpa > 0)
1235         {
1236             beta = 0;
1237             step--;   /* Don't count this step since we are restarting */
1238             continue; /* Go back to the beginning of the big for-loop */
1239         }
1240
1241         /* Calculate minimum allowed stepsize, before the average (norm)
1242          * relative change in coordinate is smaller than precision
1243          */
1244         minstep = 0;
1245         auto s_min_x = makeArrayRef(s_min->s.x);
1246         for (int i = 0; i < mdatoms->homenr; i++)
1247         {
1248             for (m = 0; m < DIM; m++)
1249             {
1250                 tmp = fabs(s_min_x[i][m]);
1251                 if (tmp < 1.0)
1252                 {
1253                     tmp = 1.0;
1254                 }
1255                 tmp      = pm[i][m]/tmp;
1256                 minstep += tmp*tmp;
1257             }
1258         }
1259         /* Add up from all CPUs */
1260         if (PAR(cr))
1261         {
1262             gmx_sumd(1, &minstep, cr);
1263         }
1264
1265         minstep = GMX_REAL_EPS/sqrt(minstep/(3*top_global->natoms));
1266
1267         if (stepsize < minstep)
1268         {
1269             converged = TRUE;
1270             break;
1271         }
1272
1273         /* Write coordinates if necessary */
1274         do_x = do_per_step(step, inputrec->nstxout);
1275         do_f = do_per_step(step, inputrec->nstfout);
1276
1277         write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
1278                       top_global, inputrec, step,
1279                       s_min, state_global, observablesHistory);
1280
1281         /* Take a step downhill.
1282          * In theory, we should minimize the function along this direction.
1283          * That is quite possible, but it turns out to take 5-10 function evaluations
1284          * for each line. However, we dont really need to find the exact minimum -
1285          * it is much better to start a new CG step in a modified direction as soon
1286          * as we are close to it. This will save a lot of energy evaluations.
1287          *
1288          * In practice, we just try to take a single step.
1289          * If it worked (i.e. lowered the energy), we increase the stepsize but
1290          * the continue straight to the next CG step without trying to find any minimum.
1291          * If it didn't work (higher energy), there must be a minimum somewhere between
1292          * the old position and the new one.
1293          *
1294          * Due to the finite numerical accuracy, it turns out that it is a good idea
1295          * to even accept a SMALL increase in energy, if the derivative is still downhill.
1296          * This leads to lower final energies in the tests I've done. / Erik
1297          */
1298         s_a->epot = s_min->epot;
1299         a         = 0.0;
1300         c         = a + stepsize; /* reference position along line is zero */
1301
1302         if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1303         {
1304             em_dd_partition_system(fplog, mdlog, step, cr, top_global, inputrec,
1305                                    s_min, &top, mdAtoms, fr, vsite, constr,
1306                                    nrnb, wcycle);
1307         }
1308
1309         /* Take a trial step (new coords in s_c) */
1310         do_em_step(cr, inputrec, mdatoms, s_min, c, &s_min->s.cg_p, s_c,
1311                    constr, -1);
1312
1313         neval++;
1314         /* Calculate energy for the trial step */
1315         energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE);
1316
1317         /* Calc derivative along line */
1318         const rvec *pc  = s_c->s.cg_p.rvec_array();
1319         const rvec *sfc = s_c->f.rvec_array();
1320         double      gpc = 0;
1321         for (int i = 0; i < mdatoms->homenr; i++)
1322         {
1323             for (m = 0; m < DIM; m++)
1324             {
1325                 gpc -= pc[i][m]*sfc[i][m]; /* f is negative gradient, thus the sign */
1326             }
1327         }
1328         /* Sum the gradient along the line across CPUs */
1329         if (PAR(cr))
1330         {
1331             gmx_sumd(1, &gpc, cr);
1332         }
1333
1334         /* This is the max amount of increase in energy we tolerate */
1335         tmp = std::sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1336
1337         /* Accept the step if the energy is lower, or if it is not significantly higher
1338          * and the line derivative is still negative.
1339          */
1340         if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1341         {
1342             foundlower = TRUE;
1343             /* Great, we found a better energy. Increase step for next iteration
1344              * if we are still going down, decrease it otherwise
1345              */
1346             if (gpc < 0)
1347             {
1348                 stepsize *= 1.618034; /* The golden section */
1349             }
1350             else
1351             {
1352                 stepsize *= 0.618034; /* 1/golden section */
1353             }
1354         }
1355         else
1356         {
1357             /* New energy is the same or higher. We will have to do some work
1358              * to find a smaller value in the interval. Take smaller step next time!
1359              */
1360             foundlower = FALSE;
1361             stepsize  *= 0.618034;
1362         }
1363
1364
1365
1366
1367         /* OK, if we didn't find a lower value we will have to locate one now - there must
1368          * be one in the interval [a=0,c].
1369          * The same thing is valid here, though: Don't spend dozens of iterations to find
1370          * the line minimum. We try to interpolate based on the derivative at the endpoints,
1371          * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1372          *
1373          * I also have a safeguard for potentially really pathological functions so we never
1374          * take more than 20 steps before we give up ...
1375          *
1376          * If we already found a lower value we just skip this step and continue to the update.
1377          */
1378         double gpb;
1379         if (!foundlower)
1380         {
1381             nminstep = 0;
1382
1383             do
1384             {
1385                 /* Select a new trial point.
1386                  * If the derivatives at points a & c have different sign we interpolate to zero,
1387                  * otherwise just do a bisection.
1388                  */
1389                 if (gpa < 0 && gpc > 0)
1390                 {
1391                     b = a + gpa*(a-c)/(gpc-gpa);
1392                 }
1393                 else
1394                 {
1395                     b = 0.5*(a+c);
1396                 }
1397
1398                 /* safeguard if interpolation close to machine accuracy causes errors:
1399                  * never go outside the interval
1400                  */
1401                 if (b <= a || b >= c)
1402                 {
1403                     b = 0.5*(a+c);
1404                 }
1405
1406                 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1407                 {
1408                     /* Reload the old state */
1409                     em_dd_partition_system(fplog, mdlog, -1, cr, top_global, inputrec,
1410                                            s_min, &top, mdAtoms, fr, vsite, constr,
1411                                            nrnb, wcycle);
1412                 }
1413
1414                 /* Take a trial step to this new point - new coords in s_b */
1415                 do_em_step(cr, inputrec, mdatoms, s_min, b, &s_min->s.cg_p, s_b,
1416                            constr, -1);
1417
1418                 neval++;
1419                 /* Calculate energy for the trial step */
1420                 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE);
1421
1422                 /* p does not change within a step, but since the domain decomposition
1423                  * might change, we have to use cg_p of s_b here.
1424                  */
1425                 const rvec *pb  = s_b->s.cg_p.rvec_array();
1426                 const rvec *sfb = s_b->f.rvec_array();
1427                 gpb             = 0;
1428                 for (int i = 0; i < mdatoms->homenr; i++)
1429                 {
1430                     for (m = 0; m < DIM; m++)
1431                     {
1432                         gpb -= pb[i][m]*sfb[i][m]; /* f is negative gradient, thus the sign */
1433                     }
1434                 }
1435                 /* Sum the gradient along the line across CPUs */
1436                 if (PAR(cr))
1437                 {
1438                     gmx_sumd(1, &gpb, cr);
1439                 }
1440
1441                 if (debug)
1442                 {
1443                     fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1444                             s_a->epot, s_b->epot, s_c->epot, gpb);
1445                 }
1446
1447                 epot_repl = s_b->epot;
1448
1449                 /* Keep one of the intervals based on the value of the derivative at the new point */
1450                 if (gpb > 0)
1451                 {
1452                     /* Replace c endpoint with b */
1453                     swap_em_state(&s_b, &s_c);
1454                     c   = b;
1455                     gpc = gpb;
1456                 }
1457                 else
1458                 {
1459                     /* Replace a endpoint with b */
1460                     swap_em_state(&s_b, &s_a);
1461                     a   = b;
1462                     gpa = gpb;
1463                 }
1464
1465                 /*
1466                  * Stop search as soon as we find a value smaller than the endpoints.
1467                  * Never run more than 20 steps, no matter what.
1468                  */
1469                 nminstep++;
1470             }
1471             while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1472                    (nminstep < 20));
1473
1474             if (std::fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1475                 nminstep >= 20)
1476             {
1477                 /* OK. We couldn't find a significantly lower energy.
1478                  * If beta==0 this was steepest descent, and then we give up.
1479                  * If not, set beta=0 and restart with steepest descent before quitting.
1480                  */
1481                 if (beta == 0.0)
1482                 {
1483                     /* Converged */
1484                     converged = TRUE;
1485                     break;
1486                 }
1487                 else
1488                 {
1489                     /* Reset memory before giving up */
1490                     beta = 0.0;
1491                     continue;
1492                 }
1493             }
1494
1495             /* Select min energy state of A & C, put the best in B.
1496              */
1497             if (s_c->epot < s_a->epot)
1498             {
1499                 if (debug)
1500                 {
1501                     fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1502                             s_c->epot, s_a->epot);
1503                 }
1504                 swap_em_state(&s_b, &s_c);
1505                 gpb = gpc;
1506             }
1507             else
1508             {
1509                 if (debug)
1510                 {
1511                     fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1512                             s_a->epot, s_c->epot);
1513                 }
1514                 swap_em_state(&s_b, &s_a);
1515                 gpb = gpa;
1516             }
1517
1518         }
1519         else
1520         {
1521             if (debug)
1522             {
1523                 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1524                         s_c->epot);
1525             }
1526             swap_em_state(&s_b, &s_c);
1527             gpb = gpc;
1528         }
1529
1530         /* new search direction */
1531         /* beta = 0 means forget all memory and restart with steepest descents. */
1532         if (nstcg && ((step % nstcg) == 0))
1533         {
1534             beta = 0.0;
1535         }
1536         else
1537         {
1538             /* s_min->fnorm cannot be zero, because then we would have converged
1539              * and broken out.
1540              */
1541
1542             /* Polak-Ribiere update.
1543              * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1544              */
1545             beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1546         }
1547         /* Limit beta to prevent oscillations */
1548         if (fabs(beta) > 5.0)
1549         {
1550             beta = 0.0;
1551         }
1552
1553
1554         /* update positions */
1555         swap_em_state(&s_min, &s_b);
1556         gpa = gpb;
1557
1558         /* Print it if necessary */
1559         if (MASTER(cr))
1560         {
1561             if (mdrunOptions.verbose)
1562             {
1563                 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1564                 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1565                         step, s_min->epot, s_min->fnorm/sqrtNumAtoms,
1566                         s_min->fmax, s_min->a_fmax+1);
1567                 fflush(stderr);
1568             }
1569             /* Store the new (lower) energies */
1570             matrix nullBox = {};
1571             energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step),
1572                                              mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1573                                              nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1574
1575             do_log = do_per_step(step, inputrec->nstlog);
1576             do_ene = do_per_step(step, inputrec->nstenergy);
1577
1578             /* Prepare IMD energy record, if bIMD is TRUE. */
1579             IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1580
1581             if (do_log)
1582             {
1583                 print_ebin_header(fplog, step, step);
1584             }
1585             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1586                                                do_log ? fplog : nullptr, step, step, eprNORMAL,
1587                                                fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1588         }
1589
1590         /* Send energies and positions to the IMD client if bIMD is TRUE. */
1591         if (MASTER(cr) && do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x.rvec_array(), inputrec, 0, wcycle))
1592         {
1593             IMD_send_positions(inputrec->imd);
1594         }
1595
1596         /* Stop when the maximum force lies below tolerance.
1597          * If we have reached machine precision, converged is already set to true.
1598          */
1599         converged = converged || (s_min->fmax < inputrec->em_tol);
1600
1601     }   /* End of the loop */
1602
1603     /* IMD cleanup, if bIMD is TRUE. */
1604     IMD_finalize(inputrec->bIMD, inputrec->imd);
1605
1606     if (converged)
1607     {
1608         step--; /* we never took that last step in this case */
1609
1610     }
1611     if (s_min->fmax > inputrec->em_tol)
1612     {
1613         if (MASTER(cr))
1614         {
1615             warn_step(fplog, inputrec->em_tol, s_min->fmax,
1616                       step-1 == number_steps, FALSE);
1617         }
1618         converged = FALSE;
1619     }
1620
1621     if (MASTER(cr))
1622     {
1623         /* If we printed energy and/or logfile last step (which was the last step)
1624          * we don't have to do it again, but otherwise print the final values.
1625          */
1626         if (!do_log)
1627         {
1628             /* Write final value to log since we didn't do anything the last step */
1629             print_ebin_header(fplog, step, step);
1630         }
1631         if (!do_ene || !do_log)
1632         {
1633             /* Write final energy file entries */
1634             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1635                                                !do_log ? fplog : nullptr, step, step, eprNORMAL,
1636                                                fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1637         }
1638     }
1639
1640     /* Print some stuff... */
1641     if (MASTER(cr))
1642     {
1643         fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1644     }
1645
1646     /* IMPORTANT!
1647      * For accurate normal mode calculation it is imperative that we
1648      * store the last conformation into the full precision binary trajectory.
1649      *
1650      * However, we should only do it if we did NOT already write this step
1651      * above (which we did if do_x or do_f was true).
1652      */
1653     /* Note that with 0 < nstfout != nstxout we can end up with two frames
1654      * in the trajectory with the same step number.
1655      */
1656     do_x = !do_per_step(step, inputrec->nstxout);
1657     do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1658
1659     write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1660                   top_global, inputrec, step,
1661                   s_min, state_global, observablesHistory);
1662
1663
1664     if (MASTER(cr))
1665     {
1666         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1667         print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1668                         s_min, sqrtNumAtoms);
1669         print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1670                         s_min, sqrtNumAtoms);
1671
1672         fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1673     }
1674
1675     finish_em(cr, outf, walltime_accounting, wcycle);
1676
1677     /* To print the actual number of steps we needed somewhere */
1678     walltime_accounting_set_nsteps_done(walltime_accounting, step);
1679 }
1680
1681
1682 void
1683 Integrator::do_lbfgs()
1684 {
1685     static const char *LBFGS = "Low-Memory BFGS Minimizer";
1686     em_state_t         ems;
1687     gmx_localtop_t     top;
1688     gmx_global_stat_t  gstat;
1689     t_graph           *graph;
1690     int                ncorr, nmaxcorr, point, cp, neval, nminstep;
1691     double             stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1692     real              *rho, *alpha, *p, *s, **dx, **dg;
1693     real               a, b, c, maxdelta, delta;
1694     real               diag, Epot0;
1695     real               dgdx, dgdg, sq, yr, beta;
1696     gmx_bool           converged;
1697     rvec               mu_tot = {0};
1698     gmx_bool           do_log, do_ene, do_x, do_f, foundlower, *frozen;
1699     tensor             vir, pres;
1700     int                start, end, number_steps;
1701     int                i, k, m, n, gf, step;
1702     int                mdof_flags;
1703     auto               mdatoms = mdAtoms->mdatoms();
1704
1705     GMX_LOG(mdlog.info).asParagraph().
1706         appendText("Note that activating L-BFGS energy minimization via the "
1707                    "integrator .mdp option and the command gmx mdrun may "
1708                    "be available in a different form in a future version of GROMACS, "
1709                    "e.g. gmx minimize and an .mdp option.");
1710
1711     if (PAR(cr))
1712     {
1713         gmx_fatal(FARGS, "L-BFGS minimization only supports a single rank");
1714     }
1715
1716     if (nullptr != constr)
1717     {
1718         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).");
1719     }
1720
1721     n        = 3*state_global->natoms;
1722     nmaxcorr = inputrec->nbfgscorr;
1723
1724     snew(frozen, n);
1725
1726     snew(p, n);
1727     snew(rho, nmaxcorr);
1728     snew(alpha, nmaxcorr);
1729
1730     snew(dx, nmaxcorr);
1731     for (i = 0; i < nmaxcorr; i++)
1732     {
1733         snew(dx[i], n);
1734     }
1735
1736     snew(dg, nmaxcorr);
1737     for (i = 0; i < nmaxcorr; i++)
1738     {
1739         snew(dg[i], n);
1740     }
1741
1742     step  = 0;
1743     neval = 0;
1744
1745     /* Init em */
1746     init_em(fplog, mdlog, LBFGS, cr, ms, inputrec, mdrunOptions,
1747             state_global, top_global, &ems, &top,
1748             nrnb, fr, &graph, mdAtoms, &gstat,
1749             vsite, constr, nullptr,
1750             nfile, fnm);
1751     gmx_mdoutf       *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, inputrec, top_global, nullptr, wcycle);
1752     gmx::EnergyOutput energyOutput;
1753     energyOutput.prepare(mdoutf_get_fp_ene(outf), top_global, inputrec, nullptr);
1754
1755     start = 0;
1756     end   = mdatoms->homenr;
1757
1758     /* We need 4 working states */
1759     em_state_t  s0 {}, s1 {}, s2 {}, s3 {};
1760     em_state_t *sa   = &s0;
1761     em_state_t *sb   = &s1;
1762     em_state_t *sc   = &s2;
1763     em_state_t *last = &s3;
1764     /* Initialize by copying the state from ems (we could skip x and f here) */
1765     *sa              = ems;
1766     *sb              = ems;
1767     *sc              = ems;
1768
1769     /* Print to log file */
1770     print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1771
1772     do_log = do_ene = do_x = do_f = TRUE;
1773
1774     /* Max number of steps */
1775     number_steps = inputrec->nsteps;
1776
1777     /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1778     gf = 0;
1779     for (i = start; i < end; i++)
1780     {
1781         if (mdatoms->cFREEZE)
1782         {
1783             gf = mdatoms->cFREEZE[i];
1784         }
1785         for (m = 0; m < DIM; m++)
1786         {
1787             frozen[3*i+m] = (inputrec->opts.nFreeze[gf][m] != 0);
1788         }
1789     }
1790     if (MASTER(cr))
1791     {
1792         sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1793     }
1794     if (fplog)
1795     {
1796         sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1797     }
1798
1799     if (vsite)
1800     {
1801         construct_vsites(vsite, state_global->x.rvec_array(), 1, nullptr,
1802                          top.idef.iparams, top.idef.il,
1803                          fr->ePBC, fr->bMolPBC, cr, state_global->box);
1804     }
1805
1806     /* Call the force routine and some auxiliary (neighboursearching etc.) */
1807     /* do_force always puts the charge groups in the box and shifts again
1808      * We do not unshift, so molecules are always whole
1809      */
1810     neval++;
1811     EnergyEvaluator energyEvaluator {
1812         fplog, mdlog, cr, ms,
1813         top_global, &top,
1814         inputrec, nrnb, wcycle, gstat,
1815         vsite, constr, fcd, graph,
1816         mdAtoms, fr, ppForceWorkload, enerd
1817     };
1818     energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE);
1819
1820     if (MASTER(cr))
1821     {
1822         /* Copy stuff to the energy bin for easy printing etc. */
1823         matrix nullBox = {};
1824         energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step),
1825                                          mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1826                                          nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1827
1828         print_ebin_header(fplog, step, step);
1829         energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE,
1830                                            fplog, step, step, eprNORMAL,
1831                                            fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1832     }
1833
1834     /* Set the initial step.
1835      * since it will be multiplied by the non-normalized search direction
1836      * vector (force vector the first time), we scale it by the
1837      * norm of the force.
1838      */
1839
1840     if (MASTER(cr))
1841     {
1842         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1843         fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1844         fprintf(stderr, "   F-max             = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1845         fprintf(stderr, "   F-Norm            = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1846         fprintf(stderr, "\n");
1847         /* and copy to the log file too... */
1848         fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1849         fprintf(fplog, "   F-max             = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1850         fprintf(fplog, "   F-Norm            = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1851         fprintf(fplog, "\n");
1852     }
1853
1854     // Point is an index to the memory of search directions, where 0 is the first one.
1855     point = 0;
1856
1857     // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1858     real *fInit = static_cast<real *>(ems.f.rvec_array()[0]);
1859     for (i = 0; i < n; i++)
1860     {
1861         if (!frozen[i])
1862         {
1863             dx[point][i] = fInit[i]; /* Initial search direction */
1864         }
1865         else
1866         {
1867             dx[point][i] = 0;
1868         }
1869     }
1870
1871     // Stepsize will be modified during the search, and actually it is not critical
1872     // (the main efficiency in the algorithm comes from changing directions), but
1873     // we still need an initial value, so estimate it as the inverse of the norm
1874     // so we take small steps where the potential fluctuates a lot.
1875     stepsize  = 1.0/ems.fnorm;
1876
1877     /* Start the loop over BFGS steps.
1878      * Each successful step is counted, and we continue until
1879      * we either converge or reach the max number of steps.
1880      */
1881
1882     ncorr = 0;
1883
1884     /* Set the gradient from the force */
1885     converged = FALSE;
1886     for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1887     {
1888
1889         /* Write coordinates if necessary */
1890         do_x = do_per_step(step, inputrec->nstxout);
1891         do_f = do_per_step(step, inputrec->nstfout);
1892
1893         mdof_flags = 0;
1894         if (do_x)
1895         {
1896             mdof_flags |= MDOF_X;
1897         }
1898
1899         if (do_f)
1900         {
1901             mdof_flags |= MDOF_F;
1902         }
1903
1904         if (inputrec->bIMD)
1905         {
1906             mdof_flags |= MDOF_IMD;
1907         }
1908
1909         mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1910                                          top_global, step, static_cast<real>(step), &ems.s, state_global, observablesHistory, ems.f);
1911
1912         /* Do the linesearching in the direction dx[point][0..(n-1)] */
1913
1914         /* make s a pointer to current search direction - point=0 first time we get here */
1915         s = dx[point];
1916
1917         real *xx = static_cast<real *>(ems.s.x.rvec_array()[0]);
1918         real *ff = static_cast<real *>(ems.f.rvec_array()[0]);
1919
1920         // calculate line gradient in position A
1921         for (gpa = 0, i = 0; i < n; i++)
1922         {
1923             gpa -= s[i]*ff[i];
1924         }
1925
1926         /* Calculate minimum allowed stepsize along the line, before the average (norm)
1927          * relative change in coordinate is smaller than precision
1928          */
1929         for (minstep = 0, i = 0; i < n; i++)
1930         {
1931             tmp = fabs(xx[i]);
1932             if (tmp < 1.0)
1933             {
1934                 tmp = 1.0;
1935             }
1936             tmp      = s[i]/tmp;
1937             minstep += tmp*tmp;
1938         }
1939         minstep = GMX_REAL_EPS/sqrt(minstep/n);
1940
1941         if (stepsize < minstep)
1942         {
1943             converged = TRUE;
1944             break;
1945         }
1946
1947         // Before taking any steps along the line, store the old position
1948         *last       = ems;
1949         real *lastx = static_cast<real *>(last->s.x.data()[0]);
1950         real *lastf = static_cast<real *>(last->f.data()[0]);
1951         Epot0       = ems.epot;
1952
1953         *sa         = ems;
1954
1955         /* Take a step downhill.
1956          * In theory, we should find the actual minimum of the function in this
1957          * direction, somewhere along the line.
1958          * That is quite possible, but it turns out to take 5-10 function evaluations
1959          * for each line. However, we dont really need to find the exact minimum -
1960          * it is much better to start a new BFGS step in a modified direction as soon
1961          * as we are close to it. This will save a lot of energy evaluations.
1962          *
1963          * In practice, we just try to take a single step.
1964          * If it worked (i.e. lowered the energy), we increase the stepsize but
1965          * continue straight to the next BFGS step without trying to find any minimum,
1966          * i.e. we change the search direction too. If the line was smooth, it is
1967          * likely we are in a smooth region, and then it makes sense to take longer
1968          * steps in the modified search direction too.
1969          *
1970          * If it didn't work (higher energy), there must be a minimum somewhere between
1971          * the old position and the new one. Then we need to start by finding a lower
1972          * value before we change search direction. Since the energy was apparently
1973          * quite rough, we need to decrease the step size.
1974          *
1975          * Due to the finite numerical accuracy, it turns out that it is a good idea
1976          * to accept a SMALL increase in energy, if the derivative is still downhill.
1977          * This leads to lower final energies in the tests I've done. / Erik
1978          */
1979
1980         // State "A" is the first position along the line.
1981         // reference position along line is initially zero
1982         a          = 0.0;
1983
1984         // Check stepsize first. We do not allow displacements
1985         // larger than emstep.
1986         //
1987         do
1988         {
1989             // Pick a new position C by adding stepsize to A.
1990             c        = a + stepsize;
1991
1992             // Calculate what the largest change in any individual coordinate
1993             // would be (translation along line * gradient along line)
1994             maxdelta = 0;
1995             for (i = 0; i < n; i++)
1996             {
1997                 delta = c*s[i];
1998                 if (delta > maxdelta)
1999                 {
2000                     maxdelta = delta;
2001                 }
2002             }
2003             // If any displacement is larger than the stepsize limit, reduce the step
2004             if (maxdelta > inputrec->em_stepsize)
2005             {
2006                 stepsize *= 0.1;
2007             }
2008         }
2009         while (maxdelta > inputrec->em_stepsize);
2010
2011         // Take a trial step and move the coordinate array xc[] to position C
2012         real *xc = static_cast<real *>(sc->s.x.rvec_array()[0]);
2013         for (i = 0; i < n; i++)
2014         {
2015             xc[i] = lastx[i] + c*s[i];
2016         }
2017
2018         neval++;
2019         // Calculate energy for the trial step in position C
2020         energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE);
2021
2022         // Calc line gradient in position C
2023         real *fc = static_cast<real *>(sc->f.rvec_array()[0]);
2024         for (gpc = 0, i = 0; i < n; i++)
2025         {
2026             gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
2027         }
2028         /* Sum the gradient along the line across CPUs */
2029         if (PAR(cr))
2030         {
2031             gmx_sumd(1, &gpc, cr);
2032         }
2033
2034         // This is the max amount of increase in energy we tolerate.
2035         // By allowing VERY small changes (close to numerical precision) we
2036         // frequently find even better (lower) final energies.
2037         tmp = std::sqrt(GMX_REAL_EPS)*fabs(sa->epot);
2038
2039         // Accept the step if the energy is lower in the new position C (compared to A),
2040         // or if it is not significantly higher and the line derivative is still negative.
2041         foundlower = sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp));
2042         // If true, great, we found a better energy. We no longer try to alter the
2043         // stepsize, but simply accept this new better position. The we select a new
2044         // search direction instead, which will be much more efficient than continuing
2045         // to take smaller steps along a line. Set fnorm based on the new C position,
2046         // which will be used to update the stepsize to 1/fnorm further down.
2047
2048         // If false, the energy is NOT lower in point C, i.e. it will be the same
2049         // or higher than in point A. In this case it is pointless to move to point C,
2050         // so we will have to do more iterations along the same line to find a smaller
2051         // value in the interval [A=0.0,C].
2052         // Here, A is still 0.0, but that will change when we do a search in the interval
2053         // [0.0,C] below. That search we will do by interpolation or bisection rather
2054         // than with the stepsize, so no need to modify it. For the next search direction
2055         // it will be reset to 1/fnorm anyway.
2056
2057         if (!foundlower)
2058         {
2059             // OK, if we didn't find a lower value we will have to locate one now - there must
2060             // be one in the interval [a,c].
2061             // The same thing is valid here, though: Don't spend dozens of iterations to find
2062             // the line minimum. We try to interpolate based on the derivative at the endpoints,
2063             // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2064             // I also have a safeguard for potentially really pathological functions so we never
2065             // take more than 20 steps before we give up.
2066             // If we already found a lower value we just skip this step and continue to the update.
2067             real fnorm = 0;
2068             nminstep   = 0;
2069             do
2070             {
2071                 // Select a new trial point B in the interval [A,C].
2072                 // If the derivatives at points a & c have different sign we interpolate to zero,
2073                 // otherwise just do a bisection since there might be multiple minima/maxima
2074                 // inside the interval.
2075                 if (gpa < 0 && gpc > 0)
2076                 {
2077                     b = a + gpa*(a-c)/(gpc-gpa);
2078                 }
2079                 else
2080                 {
2081                     b = 0.5*(a+c);
2082                 }
2083
2084                 /* safeguard if interpolation close to machine accuracy causes errors:
2085                  * never go outside the interval
2086                  */
2087                 if (b <= a || b >= c)
2088                 {
2089                     b = 0.5*(a+c);
2090                 }
2091
2092                 // Take a trial step to point B
2093                 real *xb = static_cast<real *>(sb->s.x.rvec_array()[0]);
2094                 for (i = 0; i < n; i++)
2095                 {
2096                     xb[i] = lastx[i] + b*s[i];
2097                 }
2098
2099                 neval++;
2100                 // Calculate energy for the trial step in point B
2101                 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE);
2102                 fnorm = sb->fnorm;
2103
2104                 // Calculate gradient in point B
2105                 real *fb = static_cast<real *>(sb->f.rvec_array()[0]);
2106                 for (gpb = 0, i = 0; i < n; i++)
2107                 {
2108                     gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2109
2110                 }
2111                 /* Sum the gradient along the line across CPUs */
2112                 if (PAR(cr))
2113                 {
2114                     gmx_sumd(1, &gpb, cr);
2115                 }
2116
2117                 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2118                 // at the new point B, and rename the endpoints of this new interval A and C.
2119                 if (gpb > 0)
2120                 {
2121                     /* Replace c endpoint with b */
2122                     c   = b;
2123                     /* swap states b and c */
2124                     swap_em_state(&sb, &sc);
2125                 }
2126                 else
2127                 {
2128                     /* Replace a endpoint with b */
2129                     a   = b;
2130                     /* swap states a and b */
2131                     swap_em_state(&sa, &sb);
2132                 }
2133
2134                 /*
2135                  * Stop search as soon as we find a value smaller than the endpoints,
2136                  * or if the tolerance is below machine precision.
2137                  * Never run more than 20 steps, no matter what.
2138                  */
2139                 nminstep++;
2140             }
2141             while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2142
2143             if (std::fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2144             {
2145                 /* OK. We couldn't find a significantly lower energy.
2146                  * If ncorr==0 this was steepest descent, and then we give up.
2147                  * If not, reset memory to restart as steepest descent before quitting.
2148                  */
2149                 if (ncorr == 0)
2150                 {
2151                     /* Converged */
2152                     converged = TRUE;
2153                     break;
2154                 }
2155                 else
2156                 {
2157                     /* Reset memory */
2158                     ncorr = 0;
2159                     /* Search in gradient direction */
2160                     for (i = 0; i < n; i++)
2161                     {
2162                         dx[point][i] = ff[i];
2163                     }
2164                     /* Reset stepsize */
2165                     stepsize = 1.0/fnorm;
2166                     continue;
2167                 }
2168             }
2169
2170             /* Select min energy state of A & C, put the best in xx/ff/Epot
2171              */
2172             if (sc->epot < sa->epot)
2173             {
2174                 /* Use state C */
2175                 ems        = *sc;
2176                 step_taken = c;
2177             }
2178             else
2179             {
2180                 /* Use state A */
2181                 ems        = *sa;
2182                 step_taken = a;
2183             }
2184
2185         }
2186         else
2187         {
2188             /* found lower */
2189             /* Use state C */
2190             ems        = *sc;
2191             step_taken = c;
2192         }
2193
2194         /* Update the memory information, and calculate a new
2195          * approximation of the inverse hessian
2196          */
2197
2198         /* Have new data in Epot, xx, ff */
2199         if (ncorr < nmaxcorr)
2200         {
2201             ncorr++;
2202         }
2203
2204         for (i = 0; i < n; i++)
2205         {
2206             dg[point][i]  = lastf[i]-ff[i];
2207             dx[point][i] *= step_taken;
2208         }
2209
2210         dgdg = 0;
2211         dgdx = 0;
2212         for (i = 0; i < n; i++)
2213         {
2214             dgdg += dg[point][i]*dg[point][i];
2215             dgdx += dg[point][i]*dx[point][i];
2216         }
2217
2218         diag = dgdx/dgdg;
2219
2220         rho[point] = 1.0/dgdx;
2221         point++;
2222
2223         if (point >= nmaxcorr)
2224         {
2225             point = 0;
2226         }
2227
2228         /* Update */
2229         for (i = 0; i < n; i++)
2230         {
2231             p[i] = ff[i];
2232         }
2233
2234         cp = point;
2235
2236         /* Recursive update. First go back over the memory points */
2237         for (k = 0; k < ncorr; k++)
2238         {
2239             cp--;
2240             if (cp < 0)
2241             {
2242                 cp = ncorr-1;
2243             }
2244
2245             sq = 0;
2246             for (i = 0; i < n; i++)
2247             {
2248                 sq += dx[cp][i]*p[i];
2249             }
2250
2251             alpha[cp] = rho[cp]*sq;
2252
2253             for (i = 0; i < n; i++)
2254             {
2255                 p[i] -= alpha[cp]*dg[cp][i];
2256             }
2257         }
2258
2259         for (i = 0; i < n; i++)
2260         {
2261             p[i] *= diag;
2262         }
2263
2264         /* And then go forward again */
2265         for (k = 0; k < ncorr; k++)
2266         {
2267             yr = 0;
2268             for (i = 0; i < n; i++)
2269             {
2270                 yr += p[i]*dg[cp][i];
2271             }
2272
2273             beta = rho[cp]*yr;
2274             beta = alpha[cp]-beta;
2275
2276             for (i = 0; i < n; i++)
2277             {
2278                 p[i] += beta*dx[cp][i];
2279             }
2280
2281             cp++;
2282             if (cp >= ncorr)
2283             {
2284                 cp = 0;
2285             }
2286         }
2287
2288         for (i = 0; i < n; i++)
2289         {
2290             if (!frozen[i])
2291             {
2292                 dx[point][i] = p[i];
2293             }
2294             else
2295             {
2296                 dx[point][i] = 0;
2297             }
2298         }
2299
2300         /* Print it if necessary */
2301         if (MASTER(cr))
2302         {
2303             if (mdrunOptions.verbose)
2304             {
2305                 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2306                 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2307                         step, ems.epot, ems.fnorm/sqrtNumAtoms, ems.fmax, ems.a_fmax + 1);
2308                 fflush(stderr);
2309             }
2310             /* Store the new (lower) energies */
2311             matrix nullBox = {};
2312             energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step),
2313                                              mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
2314                                              nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2315
2316             do_log = do_per_step(step, inputrec->nstlog);
2317             do_ene = do_per_step(step, inputrec->nstenergy);
2318
2319             /* Prepare IMD energy record, if bIMD is TRUE. */
2320             IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
2321
2322             if (do_log)
2323             {
2324                 print_ebin_header(fplog, step, step);
2325             }
2326             energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2327                                                do_log ? fplog : nullptr, step, step, eprNORMAL,
2328                                                fcd, &(top_global->groups), &(inputrec->opts), nullptr);
2329         }
2330
2331         /* Send x and E to IMD client, if bIMD is TRUE. */
2332         if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x.rvec_array(), inputrec, 0, wcycle) && MASTER(cr))
2333         {
2334             IMD_send_positions(inputrec->imd);
2335         }
2336
2337         // Reset stepsize in we are doing more iterations
2338         stepsize = 1.0/ems.fnorm;
2339
2340         /* Stop when the maximum force lies below tolerance.
2341          * If we have reached machine precision, converged is already set to true.
2342          */
2343         converged = converged || (ems.fmax < inputrec->em_tol);
2344
2345     }   /* End of the loop */
2346
2347     /* IMD cleanup, if bIMD is TRUE. */
2348     IMD_finalize(inputrec->bIMD, inputrec->imd);
2349
2350     if (converged)
2351     {
2352         step--; /* we never took that last step in this case */
2353
2354     }
2355     if (ems.fmax > inputrec->em_tol)
2356     {
2357         if (MASTER(cr))
2358         {
2359             warn_step(fplog, inputrec->em_tol, ems.fmax,
2360                       step-1 == number_steps, FALSE);
2361         }
2362         converged = FALSE;
2363     }
2364
2365     /* If we printed energy and/or logfile last step (which was the last step)
2366      * we don't have to do it again, but otherwise print the final values.
2367      */
2368     if (!do_log) /* Write final value to log since we didn't do anythin last step */
2369     {
2370         print_ebin_header(fplog, step, step);
2371     }
2372     if (!do_ene || !do_log) /* Write final energy file entries */
2373     {
2374         energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2375                                            !do_log ? fplog : nullptr, step, step, eprNORMAL,
2376                                            fcd, &(top_global->groups), &(inputrec->opts), nullptr);
2377     }
2378
2379     /* Print some stuff... */
2380     if (MASTER(cr))
2381     {
2382         fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2383     }
2384
2385     /* IMPORTANT!
2386      * For accurate normal mode calculation it is imperative that we
2387      * store the last conformation into the full precision binary trajectory.
2388      *
2389      * However, we should only do it if we did NOT already write this step
2390      * above (which we did if do_x or do_f was true).
2391      */
2392     do_x = !do_per_step(step, inputrec->nstxout);
2393     do_f = !do_per_step(step, inputrec->nstfout);
2394     write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2395                   top_global, inputrec, step,
2396                   &ems, state_global, observablesHistory);
2397
2398     if (MASTER(cr))
2399     {
2400         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2401         print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2402                         number_steps, &ems, sqrtNumAtoms);
2403         print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2404                         number_steps, &ems, sqrtNumAtoms);
2405
2406         fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2407     }
2408
2409     finish_em(cr, outf, walltime_accounting, wcycle);
2410
2411     /* To print the actual number of steps we needed somewhere */
2412     walltime_accounting_set_nsteps_done(walltime_accounting, step);
2413 }
2414
2415 void
2416 Integrator::do_steep()
2417 {
2418     const char       *SD  = "Steepest Descents";
2419     gmx_localtop_t    top;
2420     gmx_global_stat_t gstat;
2421     t_graph          *graph;
2422     real              stepsize;
2423     real              ustep;
2424     gmx_bool          bDone, bAbort, do_x, do_f;
2425     tensor            vir, pres;
2426     rvec              mu_tot = {0};
2427     int               nsteps;
2428     int               count          = 0;
2429     int               steps_accepted = 0;
2430     auto              mdatoms        = mdAtoms->mdatoms();
2431
2432     GMX_LOG(mdlog.info).asParagraph().
2433         appendText("Note that activating steepest-descent energy minimization via the "
2434                    "integrator .mdp option and the command gmx mdrun may "
2435                    "be available in a different form in a future version of GROMACS, "
2436                    "e.g. gmx minimize and an .mdp option.");
2437
2438     /* Create 2 states on the stack and extract pointers that we will swap */
2439     em_state_t  s0 {}, s1 {};
2440     em_state_t *s_min = &s0;
2441     em_state_t *s_try = &s1;
2442
2443     /* Init em and store the local state in s_try */
2444     init_em(fplog, mdlog, SD, cr, ms, inputrec, mdrunOptions,
2445             state_global, top_global, s_try, &top,
2446             nrnb, fr, &graph, mdAtoms, &gstat,
2447             vsite, constr, nullptr,
2448             nfile, fnm);
2449     gmx_mdoutf       *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, inputrec, top_global, nullptr, wcycle);
2450     gmx::EnergyOutput energyOutput;
2451     energyOutput.prepare(mdoutf_get_fp_ene(outf), top_global, inputrec, nullptr);
2452
2453     /* Print to log file  */
2454     print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2455
2456     /* Set variables for stepsize (in nm). This is the largest
2457      * step that we are going to make in any direction.
2458      */
2459     ustep    = inputrec->em_stepsize;
2460     stepsize = 0;
2461
2462     /* Max number of steps  */
2463     nsteps = inputrec->nsteps;
2464
2465     if (MASTER(cr))
2466     {
2467         /* Print to the screen  */
2468         sp_header(stderr, SD, inputrec->em_tol, nsteps);
2469     }
2470     if (fplog)
2471     {
2472         sp_header(fplog, SD, inputrec->em_tol, nsteps);
2473     }
2474     EnergyEvaluator energyEvaluator {
2475         fplog, mdlog, cr, ms,
2476         top_global, &top,
2477         inputrec, nrnb, wcycle, gstat,
2478         vsite, constr, fcd, graph,
2479         mdAtoms, fr, ppForceWorkload, enerd
2480     };
2481
2482     /**** HERE STARTS THE LOOP ****
2483      * count is the counter for the number of steps
2484      * bDone will be TRUE when the minimization has converged
2485      * bAbort will be TRUE when nsteps steps have been performed or when
2486      * the stepsize becomes smaller than is reasonable for machine precision
2487      */
2488     count  = 0;
2489     bDone  = FALSE;
2490     bAbort = FALSE;
2491     while (!bDone && !bAbort)
2492     {
2493         bAbort = (nsteps >= 0) && (count == nsteps);
2494
2495         /* set new coordinates, except for first step */
2496         bool validStep = true;
2497         if (count > 0)
2498         {
2499             validStep =
2500                 do_em_step(cr, inputrec, mdatoms,
2501                            s_min, stepsize, &s_min->f, s_try,
2502                            constr, count);
2503         }
2504
2505         if (validStep)
2506         {
2507             energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0);
2508         }
2509         else
2510         {
2511             // Signal constraint error during stepping with energy=inf
2512             s_try->epot = std::numeric_limits<real>::infinity();
2513         }
2514
2515         if (MASTER(cr))
2516         {
2517             print_ebin_header(fplog, count, count);
2518         }
2519
2520         if (count == 0)
2521         {
2522             s_min->epot = s_try->epot;
2523         }
2524
2525         /* Print it if necessary  */
2526         if (MASTER(cr))
2527         {
2528             if (mdrunOptions.verbose)
2529             {
2530                 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2531                         count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2532                         ( (count == 0) || (s_try->epot < s_min->epot) ) ? '\n' : '\r');
2533                 fflush(stderr);
2534             }
2535
2536             if ( (count == 0) || (s_try->epot < s_min->epot) )
2537             {
2538                 /* Store the new (lower) energies  */
2539                 matrix nullBox = {};
2540                 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(count),
2541                                                  mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
2542                                                  nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2543
2544                 /* Prepare IMD energy record, if bIMD is TRUE. */
2545                 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2546
2547                 const bool do_dr = do_per_step(steps_accepted, inputrec->nstdisreout);
2548                 const bool do_or = do_per_step(steps_accepted, inputrec->nstorireout);
2549                 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE,
2550                                                    do_dr, do_or,
2551                                                    fplog, count, count, eprNORMAL,
2552                                                    fcd, &(top_global->groups),
2553                                                    &(inputrec->opts), nullptr);
2554                 fflush(fplog);
2555             }
2556         }
2557
2558         /* Now if the new energy is smaller than the previous...
2559          * or if this is the first step!
2560          * or if we did random steps!
2561          */
2562
2563         if ( (count == 0) || (s_try->epot < s_min->epot) )
2564         {
2565             steps_accepted++;
2566
2567             /* Test whether the convergence criterion is met...  */
2568             bDone = (s_try->fmax < inputrec->em_tol);
2569
2570             /* Copy the arrays for force, positions and energy  */
2571             /* The 'Min' array always holds the coords and forces of the minimal
2572                sampled energy  */
2573             swap_em_state(&s_min, &s_try);
2574             if (count > 0)
2575             {
2576                 ustep *= 1.2;
2577             }
2578
2579             /* Write to trn, if necessary */
2580             do_x = do_per_step(steps_accepted, inputrec->nstxout);
2581             do_f = do_per_step(steps_accepted, inputrec->nstfout);
2582             write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
2583                           top_global, inputrec, count,
2584                           s_min, state_global, observablesHistory);
2585         }
2586         else
2587         {
2588             /* If energy is not smaller make the step smaller...  */
2589             ustep *= 0.5;
2590
2591             if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2592             {
2593                 /* Reload the old state */
2594                 em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec,
2595                                        s_min, &top, mdAtoms, fr, vsite, constr,
2596                                        nrnb, wcycle);
2597             }
2598         }
2599
2600         /* Determine new step  */
2601         stepsize = ustep/s_min->fmax;
2602
2603         /* Check if stepsize is too small, with 1 nm as a characteristic length */
2604 #if GMX_DOUBLE
2605         if (count == nsteps || ustep < 1e-12)
2606 #else
2607         if (count == nsteps || ustep < 1e-6)
2608 #endif
2609         {
2610             if (MASTER(cr))
2611             {
2612                 warn_step(fplog, inputrec->em_tol, s_min->fmax,
2613                           count == nsteps, constr != nullptr);
2614             }
2615             bAbort = TRUE;
2616         }
2617
2618         /* Send IMD energies and positions, if bIMD is TRUE. */
2619         if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box,
2620                    MASTER(cr) ? state_global->x.rvec_array() : nullptr,
2621                    inputrec, 0, wcycle) &&
2622             MASTER(cr))
2623         {
2624             IMD_send_positions(inputrec->imd);
2625         }
2626
2627         count++;
2628     }   /* End of the loop  */
2629
2630     /* IMD cleanup, if bIMD is TRUE. */
2631     IMD_finalize(inputrec->bIMD, inputrec->imd);
2632
2633     /* Print some data...  */
2634     if (MASTER(cr))
2635     {
2636         fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2637     }
2638     write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout != 0, ftp2fn(efSTO, nfile, fnm),
2639                   top_global, inputrec, count,
2640                   s_min, state_global, observablesHistory);
2641
2642     if (MASTER(cr))
2643     {
2644         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2645
2646         print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2647                         s_min, sqrtNumAtoms);
2648         print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2649                         s_min, sqrtNumAtoms);
2650     }
2651
2652     finish_em(cr, outf, walltime_accounting, wcycle);
2653
2654     /* To print the actual number of steps we needed somewhere */
2655     inputrec->nsteps = count;
2656
2657     walltime_accounting_set_nsteps_done(walltime_accounting, count);
2658 }
2659
2660 void
2661 Integrator::do_nm()
2662 {
2663     const char          *NM = "Normal Mode Analysis";
2664     int                  nnodes, node;
2665     gmx_localtop_t       top;
2666     gmx_global_stat_t    gstat;
2667     t_graph             *graph;
2668     tensor               vir, pres;
2669     rvec                 mu_tot = {0};
2670     rvec                *dfdx;
2671     gmx_bool             bSparse; /* use sparse matrix storage format */
2672     size_t               sz;
2673     gmx_sparsematrix_t * sparse_matrix           = nullptr;
2674     real           *     full_matrix             = nullptr;
2675
2676     /* added with respect to mdrun */
2677     int                       row, col;
2678     real                      der_range = 10.0*std::sqrt(GMX_REAL_EPS);
2679     real                      x_min;
2680     bool                      bIsMaster = MASTER(cr);
2681     auto                      mdatoms   = mdAtoms->mdatoms();
2682
2683     GMX_LOG(mdlog.info).asParagraph().
2684         appendText("Note that activating normal-mode analysis via the integrator "
2685                    ".mdp option and the command gmx mdrun may "
2686                    "be available in a different form in a future version of GROMACS, "
2687                    "e.g. gmx normal-modes.");
2688
2689     if (constr != nullptr)
2690     {
2691         gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2692     }
2693
2694     gmx_shellfc_t *shellfc;
2695
2696     em_state_t     state_work {};
2697
2698     /* Init em and store the local state in state_minimum */
2699     init_em(fplog, mdlog, NM, cr, ms, inputrec, mdrunOptions,
2700             state_global, top_global, &state_work, &top,
2701             nrnb, fr, &graph, mdAtoms, &gstat,
2702             vsite, constr, &shellfc,
2703             nfile, fnm);
2704     gmx_mdoutf            *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, inputrec, top_global, nullptr, wcycle);
2705
2706     std::vector<int>       atom_index = get_atom_index(top_global);
2707     std::vector<gmx::RVec> fneg(atom_index.size(), {0, 0, 0});
2708     snew(dfdx, atom_index.size());
2709
2710 #if !GMX_DOUBLE
2711     if (bIsMaster)
2712     {
2713         fprintf(stderr,
2714                 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2715                 "      which MIGHT not be accurate enough for normal mode analysis.\n"
2716                 "      GROMACS now uses sparse matrix storage, so the memory requirements\n"
2717                 "      are fairly modest even if you recompile in double precision.\n\n");
2718     }
2719 #endif
2720
2721     /* Check if we can/should use sparse storage format.
2722      *
2723      * Sparse format is only useful when the Hessian itself is sparse, which it
2724      * will be when we use a cutoff.
2725      * For small systems (n<1000) it is easier to always use full matrix format, though.
2726      */
2727     if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
2728     {
2729         GMX_LOG(mdlog.warning).appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2730         bSparse = FALSE;
2731     }
2732     else if (atom_index.size() < 1000)
2733     {
2734         GMX_LOG(mdlog.warning).appendTextFormatted("Small system size (N=%zu), using full Hessian format.",
2735                                                    atom_index.size());
2736         bSparse = FALSE;
2737     }
2738     else
2739     {
2740         GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2741         bSparse = TRUE;
2742     }
2743
2744     /* Number of dimensions, based on real atoms, that is not vsites or shell */
2745     sz = DIM*atom_index.size();
2746
2747     fprintf(stderr, "Allocating Hessian memory...\n\n");
2748
2749     if (bSparse)
2750     {
2751         sparse_matrix = gmx_sparsematrix_init(sz);
2752         sparse_matrix->compressed_symmetric = TRUE;
2753     }
2754     else
2755     {
2756         snew(full_matrix, sz*sz);
2757     }
2758
2759     init_nrnb(nrnb);
2760
2761
2762     /* Write start time and temperature */
2763     print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2764
2765     /* fudge nr of steps to nr of atoms */
2766     inputrec->nsteps = atom_index.size()*2;
2767
2768     if (bIsMaster)
2769     {
2770         fprintf(stderr, "starting normal mode calculation '%s'\n%" PRId64 " steps.\n\n",
2771                 *(top_global->name), inputrec->nsteps);
2772     }
2773
2774     nnodes = cr->nnodes;
2775
2776     /* Make evaluate_energy do a single node force calculation */
2777     cr->nnodes = 1;
2778     EnergyEvaluator energyEvaluator {
2779         fplog, mdlog, cr, ms,
2780         top_global, &top,
2781         inputrec, nrnb, wcycle, gstat,
2782         vsite, constr, fcd, graph,
2783         mdAtoms, fr, ppForceWorkload, enerd
2784     };
2785     energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE);
2786     cr->nnodes = nnodes;
2787
2788     /* if forces are not small, warn user */
2789     get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
2790
2791     GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
2792     if (state_work.fmax > 1.0e-3)
2793     {
2794         GMX_LOG(mdlog.warning).appendText(
2795                 "The force is probably not small enough to "
2796                 "ensure that you are at a minimum.\n"
2797                 "Be aware that negative eigenvalues may occur\n"
2798                 "when the resulting matrix is diagonalized.");
2799     }
2800
2801     /***********************************************************
2802      *
2803      *      Loop over all pairs in matrix
2804      *
2805      *      do_force called twice. Once with positive and
2806      *      once with negative displacement
2807      *
2808      ************************************************************/
2809
2810     /* Steps are divided one by one over the nodes */
2811     bool bNS          = true;
2812     auto state_work_x = makeArrayRef(state_work.s.x);
2813     auto state_work_f = makeArrayRef(state_work.f);
2814     for (unsigned int aid = cr->nodeid; aid < atom_index.size(); aid += nnodes)
2815     {
2816         size_t atom = atom_index[aid];
2817         for (size_t d = 0; d < DIM; d++)
2818         {
2819             int64_t     step        = 0;
2820             int         force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2821             double      t           = 0;
2822
2823             x_min = state_work_x[atom][d];
2824
2825             for (unsigned int dx = 0; (dx < 2); dx++)
2826             {
2827                 if (dx == 0)
2828                 {
2829                     state_work_x[atom][d] = x_min - der_range;
2830                 }
2831                 else
2832                 {
2833                     state_work_x[atom][d] = x_min + der_range;
2834                 }
2835
2836                 /* Make evaluate_energy do a single node force calculation */
2837                 cr->nnodes = 1;
2838                 if (shellfc)
2839                 {
2840                     /* Now is the time to relax the shells */
2841                     relax_shell_flexcon(fplog,
2842                                         cr,
2843                                         ms,
2844                                         mdrunOptions.verbose,
2845                                         nullptr,
2846                                         step,
2847                                         inputrec,
2848                                         bNS,
2849                                         force_flags,
2850                                         &top,
2851                                         constr,
2852                                         enerd,
2853                                         fcd,
2854                                         &state_work.s,
2855                                         state_work.f.arrayRefWithPadding(),
2856                                         vir,
2857                                         mdatoms,
2858                                         nrnb,
2859                                         wcycle,
2860                                         graph,
2861                                         &top_global->groups,
2862                                         shellfc,
2863                                         fr,
2864                                         ppForceWorkload,
2865                                         t,
2866                                         mu_tot,
2867                                         vsite,
2868                                         DDBalanceRegionHandler(nullptr));
2869                     bNS = false;
2870                     step++;
2871                 }
2872                 else
2873                 {
2874                     energyEvaluator.run(&state_work, mu_tot, vir, pres, aid*2+dx, FALSE);
2875                 }
2876
2877                 cr->nnodes = nnodes;
2878
2879                 if (dx == 0)
2880                 {
2881                     std::copy(state_work_f.begin(), state_work_f.begin()+atom_index.size(), fneg.begin());
2882                 }
2883             }
2884
2885             /* x is restored to original */
2886             state_work_x[atom][d] = x_min;
2887
2888             for (size_t j = 0; j < atom_index.size(); j++)
2889             {
2890                 for (size_t k = 0; (k < DIM); k++)
2891                 {
2892                     dfdx[j][k] =
2893                         -(state_work_f[atom_index[j]][k] - fneg[j][k])/(2*der_range);
2894                 }
2895             }
2896
2897             if (!bIsMaster)
2898             {
2899 #if GMX_MPI
2900 #define mpi_type GMX_MPI_REAL
2901                 MPI_Send(dfdx[0], atom_index.size()*DIM, mpi_type, MASTER(cr),
2902                          cr->nodeid, cr->mpi_comm_mygroup);
2903 #endif
2904             }
2905             else
2906             {
2907                 for (node = 0; (node < nnodes && aid+node < atom_index.size()); node++)
2908                 {
2909                     if (node > 0)
2910                     {
2911 #if GMX_MPI
2912                         MPI_Status stat;
2913                         MPI_Recv(dfdx[0], atom_index.size()*DIM, mpi_type, node, node,
2914                                  cr->mpi_comm_mygroup, &stat);
2915 #undef mpi_type
2916 #endif
2917                     }
2918
2919                     row = (aid + node)*DIM + d;
2920
2921                     for (size_t j = 0; j < atom_index.size(); j++)
2922                     {
2923                         for (size_t k = 0; k < DIM; k++)
2924                         {
2925                             col = j*DIM + k;
2926
2927                             if (bSparse)
2928                             {
2929                                 if (col >= row && dfdx[j][k] != 0.0)
2930                                 {
2931                                     gmx_sparsematrix_increment_value(sparse_matrix,
2932                                                                      row, col, dfdx[j][k]);
2933                                 }
2934                             }
2935                             else
2936                             {
2937                                 full_matrix[row*sz+col] = dfdx[j][k];
2938                             }
2939                         }
2940                     }
2941                 }
2942             }
2943
2944             if (mdrunOptions.verbose && fplog)
2945             {
2946                 fflush(fplog);
2947             }
2948         }
2949         /* write progress */
2950         if (bIsMaster && mdrunOptions.verbose)
2951         {
2952             fprintf(stderr, "\rFinished step %d out of %td",
2953                     std::min<int>(atom+nnodes, atom_index.size()),
2954                     ssize(atom_index));
2955             fflush(stderr);
2956         }
2957     }
2958
2959     if (bIsMaster)
2960     {
2961         fprintf(stderr, "\n\nWriting Hessian...\n");
2962         gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2963     }
2964
2965     finish_em(cr, outf, walltime_accounting, wcycle);
2966
2967     walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size()*2);
2968 }
2969
2970 } // namespace gmx