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