Enable more warnings for Clang 6
[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         // We should move this check to the different minimizers
700         if (!validStep && ir->eI != eiSteep)
701         {
702             gmx_fatal(FARGS, "The coordinates could not be constrained. Minimizer '%s' can not handle constraint failures, use minimizer '%s' before using '%s'.",
703                       EI(ir->eI), EI(eiSteep), EI(ir->eI));
704         }
705     }
706
707     return validStep;
708 }
709
710 //! Prepare EM for using domain decomposition parallellization
711 static void em_dd_partition_system(FILE *fplog, int step, const t_commrec *cr,
712                                    gmx_mtop_t *top_global, t_inputrec *ir,
713                                    em_state_t *ems, gmx_localtop_t *top,
714                                    gmx::MDAtoms *mdAtoms, t_forcerec *fr,
715                                    gmx_vsite_t *vsite, gmx::Constraints *constr,
716                                    t_nrnb *nrnb, gmx_wallcycle_t wcycle)
717 {
718     /* Repartition the domain decomposition */
719     dd_partition_system(fplog, step, cr, FALSE, 1,
720                         nullptr, top_global, ir,
721                         &ems->s, &ems->f,
722                         mdAtoms, top, fr, vsite, constr,
723                         nrnb, wcycle, FALSE);
724     dd_store_state(cr->dd, &ems->s);
725 }
726
727 namespace
728 {
729
730 /*! \brief Class to handle the work of setting and doing an energy evaluation.
731  *
732  * This class is a mere aggregate of parameters to pass to evaluate an
733  * energy, so that future changes to names and types of them consume
734  * less time when refactoring other code.
735  *
736  * Aggregate initialization is used, for which the chief risk is that
737  * if a member is added at the end and not all initializer lists are
738  * updated, then the member will be value initialized, which will
739  * typically mean initialization to zero.
740  *
741  * We only want to construct one of these with an initializer list, so
742  * we explicitly delete the default constructor. */
743 class EnergyEvaluator
744 {
745     public:
746         //! We only intend to construct such objects with an initializer list.
747 #if __GNUC__ > 4 || (__GNUC__ == 4 && __GNUC_MINOR__ >= 9)
748         // Aspects of the C++11 spec changed after GCC 4.8.5, and
749         // compilation of the initializer list construction in
750         // runner.cpp fails in GCC 4.8.5.
751         EnergyEvaluator() = delete;
752 #endif
753         /*! \brief Evaluates an energy on the state in \c ems.
754          *
755          * \todo In practice, the same objects mu_tot, vir, and pres
756          * are always passed to this function, so we would rather have
757          * them as data members. However, their C-array types are
758          * unsuited for aggregate initialization. When the types
759          * improve, the call signature of this method can be reduced.
760          */
761         void run(em_state_t *ems, rvec mu_tot,
762                  tensor vir, tensor pres,
763                  gmx_int64_t count, gmx_bool bFirst);
764         //! Handles logging.
765         FILE                 *fplog;
766         //! Handles communication.
767         const t_commrec      *cr;
768         //! Coordinates multi-simulations.
769         const gmx_multisim_t *ms;
770         //! Holds the simulation topology.
771         gmx_mtop_t           *top_global;
772         //! Holds the domain topology.
773         gmx_localtop_t       *top;
774         //! User input options.
775         t_inputrec           *inputrec;
776         //! Manages flop accounting.
777         t_nrnb               *nrnb;
778         //! Manages wall cycle accounting.
779         gmx_wallcycle_t       wcycle;
780         //! Coordinates global reduction.
781         gmx_global_stat_t     gstat;
782         //! Handles virtual sites.
783         gmx_vsite_t          *vsite;
784         //! Handles constraints.
785         gmx::Constraints     *constr;
786         //! Handles strange things.
787         t_fcdata             *fcd;
788         //! Molecular graph for SHAKE.
789         t_graph              *graph;
790         //! Per-atom data for this domain.
791         gmx::MDAtoms         *mdAtoms;
792         //! Handles how to calculate the forces.
793         t_forcerec           *fr;
794         //! Stores the computed energies.
795         gmx_enerdata_t       *enerd;
796 };
797
798 void
799 EnergyEvaluator::run(em_state_t *ems, rvec mu_tot,
800                      tensor vir, tensor pres,
801                      gmx_int64_t count, gmx_bool bFirst)
802 {
803     real     t;
804     gmx_bool bNS;
805     tensor   force_vir, shake_vir, ekin;
806     real     dvdl_constr, prescorr, enercorr, dvdlcorr;
807     real     terminate = 0;
808
809     /* Set the time to the initial time, the time does not change during EM */
810     t = inputrec->init_t;
811
812     if (bFirst ||
813         (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
814     {
815         /* This is the first state or an old state used before the last ns */
816         bNS = TRUE;
817     }
818     else
819     {
820         bNS = FALSE;
821         if (inputrec->nstlist > 0)
822         {
823             bNS = TRUE;
824         }
825     }
826
827     if (vsite)
828     {
829         construct_vsites(vsite, as_rvec_array(ems->s.x.data()), 1, nullptr,
830                          top->idef.iparams, top->idef.il,
831                          fr->ePBC, fr->bMolPBC, cr, ems->s.box);
832     }
833
834     if (DOMAINDECOMP(cr) && bNS)
835     {
836         /* Repartition the domain decomposition */
837         em_dd_partition_system(fplog, count, cr, top_global, inputrec,
838                                ems, top, mdAtoms, fr, vsite, constr,
839                                nrnb, wcycle);
840     }
841
842     /* Calc force & energy on new trial position  */
843     /* do_force always puts the charge groups in the box and shifts again
844      * We do not unshift, so molecules are always whole in congrad.c
845      */
846     do_force(fplog, cr, ms, inputrec, nullptr,
847              count, nrnb, wcycle, top, &top_global->groups,
848              ems->s.box, ems->s.x, &ems->s.hist,
849              ems->f, force_vir, mdAtoms->mdatoms(), enerd, fcd,
850              ems->s.lambda, graph, fr, vsite, mu_tot, t, nullptr,
851              GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
852              GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
853              (bNS ? GMX_FORCE_NS : 0),
854              DOMAINDECOMP(cr) ?
855              DdOpenBalanceRegionBeforeForceComputation::yes :
856              DdOpenBalanceRegionBeforeForceComputation::no,
857              DOMAINDECOMP(cr) ?
858              DdCloseBalanceRegionAfterForceComputation::yes :
859              DdCloseBalanceRegionAfterForceComputation::no);
860
861     /* Clear the unused shake virial and pressure */
862     clear_mat(shake_vir);
863     clear_mat(pres);
864
865     /* Communicate stuff when parallel */
866     if (PAR(cr) && inputrec->eI != eiNM)
867     {
868         wallcycle_start(wcycle, ewcMoveE);
869
870         global_stat(gstat, cr, enerd, force_vir, shake_vir, mu_tot,
871                     inputrec, nullptr, nullptr, nullptr, 1, &terminate,
872                     nullptr, FALSE,
873                     CGLO_ENERGY |
874                     CGLO_PRESSURE |
875                     CGLO_CONSTRAINT);
876
877         wallcycle_stop(wcycle, ewcMoveE);
878     }
879
880     /* Calculate long range corrections to pressure and energy */
881     calc_dispcorr(inputrec, fr, ems->s.box, ems->s.lambda[efptVDW],
882                   pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
883     enerd->term[F_DISPCORR] = enercorr;
884     enerd->term[F_EPOT]    += enercorr;
885     enerd->term[F_PRES]    += prescorr;
886     enerd->term[F_DVDL]    += dvdlcorr;
887
888     ems->epot = enerd->term[F_EPOT];
889
890     if (constr)
891     {
892         /* Project out the constraint components of the force */
893         dvdl_constr = 0;
894         rvec *f_rvec = as_rvec_array(ems->f.data());
895         constr->apply(FALSE, FALSE,
896                       count, 0, 1.0,
897                       as_rvec_array(ems->s.x.data()), f_rvec, f_rvec,
898                       ems->s.box,
899                       ems->s.lambda[efptBONDED], &dvdl_constr,
900                       nullptr, &shake_vir, gmx::ConstraintVariable::ForceDispl);
901         enerd->term[F_DVDL_CONSTR] += dvdl_constr;
902         m_add(force_vir, shake_vir, vir);
903     }
904     else
905     {
906         copy_mat(force_vir, vir);
907     }
908
909     clear_mat(ekin);
910     enerd->term[F_PRES] =
911         calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
912
913     sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
914
915     if (EI_ENERGY_MINIMIZATION(inputrec->eI))
916     {
917         get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
918     }
919 }
920
921 } // namespace
922
923 //! Parallel utility summing energies and forces
924 static double reorder_partsum(const t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
925                               gmx_mtop_t *top_global,
926                               em_state_t *s_min, em_state_t *s_b)
927 {
928     t_block       *cgs_gl;
929     int            ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
930     double         partsum;
931     unsigned char *grpnrFREEZE;
932
933     if (debug)
934     {
935         fprintf(debug, "Doing reorder_partsum\n");
936     }
937
938     const rvec *fm = as_rvec_array(s_min->f.data());
939     const rvec *fb = as_rvec_array(s_b->f.data());
940
941     cgs_gl = dd_charge_groups_global(cr->dd);
942     index  = cgs_gl->index;
943
944     /* Collect fm in a global vector fmg.
945      * This conflicts with the spirit of domain decomposition,
946      * but to fully optimize this a much more complicated algorithm is required.
947      */
948     rvec *fmg;
949     snew(fmg, top_global->natoms);
950
951     ncg   = s_min->s.cg_gl.size();
952     cg_gl = s_min->s.cg_gl.data();
953     i     = 0;
954     for (c = 0; c < ncg; c++)
955     {
956         cg = cg_gl[c];
957         a0 = index[cg];
958         a1 = index[cg+1];
959         for (a = a0; a < a1; a++)
960         {
961             copy_rvec(fm[i], fmg[a]);
962             i++;
963         }
964     }
965     gmx_sum(top_global->natoms*3, fmg[0], cr);
966
967     /* Now we will determine the part of the sum for the cgs in state s_b */
968     ncg         = s_b->s.cg_gl.size();
969     cg_gl       = s_b->s.cg_gl.data();
970     partsum     = 0;
971     i           = 0;
972     gf          = 0;
973     grpnrFREEZE = top_global->groups.grpnr[egcFREEZE];
974     for (c = 0; c < ncg; c++)
975     {
976         cg = cg_gl[c];
977         a0 = index[cg];
978         a1 = index[cg+1];
979         for (a = a0; a < a1; a++)
980         {
981             if (mdatoms->cFREEZE && grpnrFREEZE)
982             {
983                 gf = grpnrFREEZE[i];
984             }
985             for (m = 0; m < DIM; m++)
986             {
987                 if (!opts->nFreeze[gf][m])
988                 {
989                     partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
990                 }
991             }
992             i++;
993         }
994     }
995
996     sfree(fmg);
997
998     return partsum;
999 }
1000
1001 //! Print some stuff, like beta, whatever that means.
1002 static real pr_beta(const t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
1003                     gmx_mtop_t *top_global,
1004                     em_state_t *s_min, em_state_t *s_b)
1005 {
1006     double sum;
1007
1008     /* This is just the classical Polak-Ribiere calculation of beta;
1009      * it looks a bit complicated since we take freeze groups into account,
1010      * and might have to sum it in parallel runs.
1011      */
1012
1013     if (!DOMAINDECOMP(cr) ||
1014         (s_min->s.ddp_count == cr->dd->ddp_count &&
1015          s_b->s.ddp_count   == cr->dd->ddp_count))
1016     {
1017         const rvec *fm  = as_rvec_array(s_min->f.data());
1018         const rvec *fb  = as_rvec_array(s_b->f.data());
1019         sum             = 0;
1020         int         gf  = 0;
1021         /* This part of code can be incorrect with DD,
1022          * since the atom ordering in s_b and s_min might differ.
1023          */
1024         for (int i = 0; i < mdatoms->homenr; i++)
1025         {
1026             if (mdatoms->cFREEZE)
1027             {
1028                 gf = mdatoms->cFREEZE[i];
1029             }
1030             for (int m = 0; m < DIM; m++)
1031             {
1032                 if (!opts->nFreeze[gf][m])
1033                 {
1034                     sum += (fb[i][m] - fm[i][m])*fb[i][m];
1035                 }
1036             }
1037         }
1038     }
1039     else
1040     {
1041         /* We need to reorder cgs while summing */
1042         sum = reorder_partsum(cr, opts, mdatoms, top_global, s_min, s_b);
1043     }
1044     if (PAR(cr))
1045     {
1046         gmx_sumd(1, &sum, cr);
1047     }
1048
1049     return sum/gmx::square(s_min->fnorm);
1050 }
1051
1052 namespace gmx
1053 {
1054
1055 void
1056 Integrator::do_cg()
1057 {
1058     const char       *CG = "Polak-Ribiere Conjugate Gradients";
1059
1060     gmx_localtop_t   *top;
1061     gmx_enerdata_t   *enerd;
1062     gmx_global_stat_t gstat;
1063     t_graph          *graph;
1064     double            tmp, minstep;
1065     real              stepsize;
1066     real              a, b, c, beta = 0.0;
1067     real              epot_repl = 0;
1068     real              pnorm;
1069     t_mdebin         *mdebin;
1070     gmx_bool          converged, foundlower;
1071     rvec              mu_tot;
1072     gmx_bool          do_log = FALSE, do_ene = FALSE, do_x, do_f;
1073     tensor            vir, pres;
1074     int               number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1075     gmx_mdoutf_t      outf;
1076     int               m, step, nminstep;
1077     auto              mdatoms = mdAtoms->mdatoms();
1078
1079     step = 0;
1080
1081     // Ensure the extra per-atom state array gets allocated
1082     state_global->flags |= (1<<estCGP);
1083
1084     /* Create 4 states on the stack and extract pointers that we will swap */
1085     em_state_t  s0 {}, s1 {}, s2 {}, s3 {};
1086     em_state_t *s_min = &s0;
1087     em_state_t *s_a   = &s1;
1088     em_state_t *s_b   = &s2;
1089     em_state_t *s_c   = &s3;
1090
1091     /* Init em and store the local state in s_min */
1092     init_em(fplog, CG, cr, ms, outputProvider, inputrec, mdrunOptions,
1093             state_global, top_global, s_min, &top,
1094             nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
1095             vsite, constr, nullptr,
1096             nfile, fnm, &outf, &mdebin, wcycle);
1097
1098     /* Print to log file */
1099     print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1100
1101     /* Max number of steps */
1102     number_steps = inputrec->nsteps;
1103
1104     if (MASTER(cr))
1105     {
1106         sp_header(stderr, CG, inputrec->em_tol, number_steps);
1107     }
1108     if (fplog)
1109     {
1110         sp_header(fplog, CG, inputrec->em_tol, number_steps);
1111     }
1112
1113     EnergyEvaluator energyEvaluator {
1114         fplog, cr, ms,
1115         top_global, top,
1116         inputrec, nrnb, wcycle, gstat,
1117         vsite, constr, fcd, graph,
1118         mdAtoms, fr, enerd
1119     };
1120     /* Call the force routine and some auxiliary (neighboursearching etc.) */
1121     /* do_force always puts the charge groups in the box and shifts again
1122      * We do not unshift, so molecules are always whole in congrad.c
1123      */
1124     energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE);
1125
1126     if (MASTER(cr))
1127     {
1128         /* Copy stuff to the energy bin for easy printing etc. */
1129         upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1130                    mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1131                    nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1132
1133         print_ebin_header(fplog, step, step);
1134         print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1135                    mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1136     }
1137
1138     /* Estimate/guess the initial stepsize */
1139     stepsize = inputrec->em_stepsize/s_min->fnorm;
1140
1141     if (MASTER(cr))
1142     {
1143         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1144         fprintf(stderr, "   F-max             = %12.5e on atom %d\n",
1145                 s_min->fmax, s_min->a_fmax+1);
1146         fprintf(stderr, "   F-Norm            = %12.5e\n",
1147                 s_min->fnorm/sqrtNumAtoms);
1148         fprintf(stderr, "\n");
1149         /* and copy to the log file too... */
1150         fprintf(fplog, "   F-max             = %12.5e on atom %d\n",
1151                 s_min->fmax, s_min->a_fmax+1);
1152         fprintf(fplog, "   F-Norm            = %12.5e\n",
1153                 s_min->fnorm/sqrtNumAtoms);
1154         fprintf(fplog, "\n");
1155     }
1156     /* Start the loop over CG steps.
1157      * Each successful step is counted, and we continue until
1158      * we either converge or reach the max number of steps.
1159      */
1160     converged = FALSE;
1161     for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1162     {
1163
1164         /* start taking steps in a new direction
1165          * First time we enter the routine, beta=0, and the direction is
1166          * simply the negative gradient.
1167          */
1168
1169         /* Calculate the new direction in p, and the gradient in this direction, gpa */
1170         rvec       *pm  = as_rvec_array(s_min->s.cg_p.data());
1171         const rvec *sfm = as_rvec_array(s_min->f.data());
1172         double      gpa = 0;
1173         int         gf  = 0;
1174         for (int i = 0; i < mdatoms->homenr; i++)
1175         {
1176             if (mdatoms->cFREEZE)
1177             {
1178                 gf = mdatoms->cFREEZE[i];
1179             }
1180             for (m = 0; m < DIM; m++)
1181             {
1182                 if (!inputrec->opts.nFreeze[gf][m])
1183                 {
1184                     pm[i][m] = sfm[i][m] + beta*pm[i][m];
1185                     gpa     -= pm[i][m]*sfm[i][m];
1186                     /* f is negative gradient, thus the sign */
1187                 }
1188                 else
1189                 {
1190                     pm[i][m] = 0;
1191                 }
1192             }
1193         }
1194
1195         /* Sum the gradient along the line across CPUs */
1196         if (PAR(cr))
1197         {
1198             gmx_sumd(1, &gpa, cr);
1199         }
1200
1201         /* Calculate the norm of the search vector */
1202         get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1203
1204         /* Just in case stepsize reaches zero due to numerical precision... */
1205         if (stepsize <= 0)
1206         {
1207             stepsize = inputrec->em_stepsize/pnorm;
1208         }
1209
1210         /*
1211          * Double check the value of the derivative in the search direction.
1212          * If it is positive it must be due to the old information in the
1213          * CG formula, so just remove that and start over with beta=0.
1214          * This corresponds to a steepest descent step.
1215          */
1216         if (gpa > 0)
1217         {
1218             beta = 0;
1219             step--;   /* Don't count this step since we are restarting */
1220             continue; /* Go back to the beginning of the big for-loop */
1221         }
1222
1223         /* Calculate minimum allowed stepsize, before the average (norm)
1224          * relative change in coordinate is smaller than precision
1225          */
1226         minstep = 0;
1227         for (int i = 0; i < mdatoms->homenr; i++)
1228         {
1229             for (m = 0; m < DIM; m++)
1230             {
1231                 tmp = fabs(s_min->s.x[i][m]);
1232                 if (tmp < 1.0)
1233                 {
1234                     tmp = 1.0;
1235                 }
1236                 tmp      = pm[i][m]/tmp;
1237                 minstep += tmp*tmp;
1238             }
1239         }
1240         /* Add up from all CPUs */
1241         if (PAR(cr))
1242         {
1243             gmx_sumd(1, &minstep, cr);
1244         }
1245
1246         minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1247
1248         if (stepsize < minstep)
1249         {
1250             converged = TRUE;
1251             break;
1252         }
1253
1254         /* Write coordinates if necessary */
1255         do_x = do_per_step(step, inputrec->nstxout);
1256         do_f = do_per_step(step, inputrec->nstfout);
1257
1258         write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
1259                       top_global, inputrec, step,
1260                       s_min, state_global, observablesHistory);
1261
1262         /* Take a step downhill.
1263          * In theory, we should minimize the function along this direction.
1264          * That is quite possible, but it turns out to take 5-10 function evaluations
1265          * for each line. However, we dont really need to find the exact minimum -
1266          * it is much better to start a new CG step in a modified direction as soon
1267          * as we are close to it. This will save a lot of energy evaluations.
1268          *
1269          * In practice, we just try to take a single step.
1270          * If it worked (i.e. lowered the energy), we increase the stepsize but
1271          * the continue straight to the next CG step without trying to find any minimum.
1272          * If it didn't work (higher energy), there must be a minimum somewhere between
1273          * the old position and the new one.
1274          *
1275          * Due to the finite numerical accuracy, it turns out that it is a good idea
1276          * to even accept a SMALL increase in energy, if the derivative is still downhill.
1277          * This leads to lower final energies in the tests I've done. / Erik
1278          */
1279         s_a->epot = s_min->epot;
1280         a         = 0.0;
1281         c         = a + stepsize; /* reference position along line is zero */
1282
1283         if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1284         {
1285             em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1286                                    s_min, top, mdAtoms, fr, vsite, constr,
1287                                    nrnb, wcycle);
1288         }
1289
1290         /* Take a trial step (new coords in s_c) */
1291         do_em_step(cr, inputrec, mdatoms, s_min, c, &s_min->s.cg_p, s_c,
1292                    constr, -1);
1293
1294         neval++;
1295         /* Calculate energy for the trial step */
1296         energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE);
1297
1298         /* Calc derivative along line */
1299         const rvec *pc  = as_rvec_array(s_c->s.cg_p.data());
1300         const rvec *sfc = as_rvec_array(s_c->f.data());
1301         double      gpc = 0;
1302         for (int i = 0; i < mdatoms->homenr; i++)
1303         {
1304             for (m = 0; m < DIM; m++)
1305             {
1306                 gpc -= pc[i][m]*sfc[i][m]; /* f is negative gradient, thus the sign */
1307             }
1308         }
1309         /* Sum the gradient along the line across CPUs */
1310         if (PAR(cr))
1311         {
1312             gmx_sumd(1, &gpc, cr);
1313         }
1314
1315         /* This is the max amount of increase in energy we tolerate */
1316         tmp = std::sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1317
1318         /* Accept the step if the energy is lower, or if it is not significantly higher
1319          * and the line derivative is still negative.
1320          */
1321         if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1322         {
1323             foundlower = TRUE;
1324             /* Great, we found a better energy. Increase step for next iteration
1325              * if we are still going down, decrease it otherwise
1326              */
1327             if (gpc < 0)
1328             {
1329                 stepsize *= 1.618034; /* The golden section */
1330             }
1331             else
1332             {
1333                 stepsize *= 0.618034; /* 1/golden section */
1334             }
1335         }
1336         else
1337         {
1338             /* New energy is the same or higher. We will have to do some work
1339              * to find a smaller value in the interval. Take smaller step next time!
1340              */
1341             foundlower = FALSE;
1342             stepsize  *= 0.618034;
1343         }
1344
1345
1346
1347
1348         /* OK, if we didn't find a lower value we will have to locate one now - there must
1349          * be one in the interval [a=0,c].
1350          * The same thing is valid here, though: Don't spend dozens of iterations to find
1351          * the line minimum. We try to interpolate based on the derivative at the endpoints,
1352          * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1353          *
1354          * I also have a safeguard for potentially really pathological functions so we never
1355          * take more than 20 steps before we give up ...
1356          *
1357          * If we already found a lower value we just skip this step and continue to the update.
1358          */
1359         double gpb;
1360         if (!foundlower)
1361         {
1362             nminstep = 0;
1363
1364             do
1365             {
1366                 /* Select a new trial point.
1367                  * If the derivatives at points a & c have different sign we interpolate to zero,
1368                  * otherwise just do a bisection.
1369                  */
1370                 if (gpa < 0 && gpc > 0)
1371                 {
1372                     b = a + gpa*(a-c)/(gpc-gpa);
1373                 }
1374                 else
1375                 {
1376                     b = 0.5*(a+c);
1377                 }
1378
1379                 /* safeguard if interpolation close to machine accuracy causes errors:
1380                  * never go outside the interval
1381                  */
1382                 if (b <= a || b >= c)
1383                 {
1384                     b = 0.5*(a+c);
1385                 }
1386
1387                 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1388                 {
1389                     /* Reload the old state */
1390                     em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1391                                            s_min, top, mdAtoms, fr, vsite, constr,
1392                                            nrnb, wcycle);
1393                 }
1394
1395                 /* Take a trial step to this new point - new coords in s_b */
1396                 do_em_step(cr, inputrec, mdatoms, s_min, b, &s_min->s.cg_p, s_b,
1397                            constr, -1);
1398
1399                 neval++;
1400                 /* Calculate energy for the trial step */
1401                 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE);
1402
1403                 /* p does not change within a step, but since the domain decomposition
1404                  * might change, we have to use cg_p of s_b here.
1405                  */
1406                 const rvec *pb  = as_rvec_array(s_b->s.cg_p.data());
1407                 const rvec *sfb = as_rvec_array(s_b->f.data());
1408                 gpb             = 0;
1409                 for (int i = 0; i < mdatoms->homenr; i++)
1410                 {
1411                     for (m = 0; m < DIM; m++)
1412                     {
1413                         gpb -= pb[i][m]*sfb[i][m]; /* f is negative gradient, thus the sign */
1414                     }
1415                 }
1416                 /* Sum the gradient along the line across CPUs */
1417                 if (PAR(cr))
1418                 {
1419                     gmx_sumd(1, &gpb, cr);
1420                 }
1421
1422                 if (debug)
1423                 {
1424                     fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1425                             s_a->epot, s_b->epot, s_c->epot, gpb);
1426                 }
1427
1428                 epot_repl = s_b->epot;
1429
1430                 /* Keep one of the intervals based on the value of the derivative at the new point */
1431                 if (gpb > 0)
1432                 {
1433                     /* Replace c endpoint with b */
1434                     swap_em_state(&s_b, &s_c);
1435                     c   = b;
1436                     gpc = gpb;
1437                 }
1438                 else
1439                 {
1440                     /* Replace a endpoint with b */
1441                     swap_em_state(&s_b, &s_a);
1442                     a   = b;
1443                     gpa = gpb;
1444                 }
1445
1446                 /*
1447                  * Stop search as soon as we find a value smaller than the endpoints.
1448                  * Never run more than 20 steps, no matter what.
1449                  */
1450                 nminstep++;
1451             }
1452             while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1453                    (nminstep < 20));
1454
1455             if (std::fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1456                 nminstep >= 20)
1457             {
1458                 /* OK. We couldn't find a significantly lower energy.
1459                  * If beta==0 this was steepest descent, and then we give up.
1460                  * If not, set beta=0 and restart with steepest descent before quitting.
1461                  */
1462                 if (beta == 0.0)
1463                 {
1464                     /* Converged */
1465                     converged = TRUE;
1466                     break;
1467                 }
1468                 else
1469                 {
1470                     /* Reset memory before giving up */
1471                     beta = 0.0;
1472                     continue;
1473                 }
1474             }
1475
1476             /* Select min energy state of A & C, put the best in B.
1477              */
1478             if (s_c->epot < s_a->epot)
1479             {
1480                 if (debug)
1481                 {
1482                     fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1483                             s_c->epot, s_a->epot);
1484                 }
1485                 swap_em_state(&s_b, &s_c);
1486                 gpb = gpc;
1487             }
1488             else
1489             {
1490                 if (debug)
1491                 {
1492                     fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1493                             s_a->epot, s_c->epot);
1494                 }
1495                 swap_em_state(&s_b, &s_a);
1496                 gpb = gpa;
1497             }
1498
1499         }
1500         else
1501         {
1502             if (debug)
1503             {
1504                 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1505                         s_c->epot);
1506             }
1507             swap_em_state(&s_b, &s_c);
1508             gpb = gpc;
1509         }
1510
1511         /* new search direction */
1512         /* beta = 0 means forget all memory and restart with steepest descents. */
1513         if (nstcg && ((step % nstcg) == 0))
1514         {
1515             beta = 0.0;
1516         }
1517         else
1518         {
1519             /* s_min->fnorm cannot be zero, because then we would have converged
1520              * and broken out.
1521              */
1522
1523             /* Polak-Ribiere update.
1524              * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1525              */
1526             beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1527         }
1528         /* Limit beta to prevent oscillations */
1529         if (fabs(beta) > 5.0)
1530         {
1531             beta = 0.0;
1532         }
1533
1534
1535         /* update positions */
1536         swap_em_state(&s_min, &s_b);
1537         gpa = gpb;
1538
1539         /* Print it if necessary */
1540         if (MASTER(cr))
1541         {
1542             if (mdrunOptions.verbose)
1543             {
1544                 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1545                 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1546                         step, s_min->epot, s_min->fnorm/sqrtNumAtoms,
1547                         s_min->fmax, s_min->a_fmax+1);
1548                 fflush(stderr);
1549             }
1550             /* Store the new (lower) energies */
1551             upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1552                        mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1553                        nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1554
1555             do_log = do_per_step(step, inputrec->nstlog);
1556             do_ene = do_per_step(step, inputrec->nstenergy);
1557
1558             /* Prepare IMD energy record, if bIMD is TRUE. */
1559             IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1560
1561             if (do_log)
1562             {
1563                 print_ebin_header(fplog, step, step);
1564             }
1565             print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1566                        do_log ? fplog : nullptr, step, step, eprNORMAL,
1567                        mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1568         }
1569
1570         /* Send energies and positions to the IMD client if bIMD is TRUE. */
1571         if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, as_rvec_array(state_global->x.data()), inputrec, 0, wcycle) && MASTER(cr))
1572         {
1573             IMD_send_positions(inputrec->imd);
1574         }
1575
1576         /* Stop when the maximum force lies below tolerance.
1577          * If we have reached machine precision, converged is already set to true.
1578          */
1579         converged = converged || (s_min->fmax < inputrec->em_tol);
1580
1581     }   /* End of the loop */
1582
1583     /* IMD cleanup, if bIMD is TRUE. */
1584     IMD_finalize(inputrec->bIMD, inputrec->imd);
1585
1586     if (converged)
1587     {
1588         step--; /* we never took that last step in this case */
1589
1590     }
1591     if (s_min->fmax > inputrec->em_tol)
1592     {
1593         if (MASTER(cr))
1594         {
1595             warn_step(fplog, inputrec->em_tol, s_min->fmax,
1596                       step-1 == number_steps, FALSE);
1597         }
1598         converged = FALSE;
1599     }
1600
1601     if (MASTER(cr))
1602     {
1603         /* If we printed energy and/or logfile last step (which was the last step)
1604          * we don't have to do it again, but otherwise print the final values.
1605          */
1606         if (!do_log)
1607         {
1608             /* Write final value to log since we didn't do anything the last step */
1609             print_ebin_header(fplog, step, step);
1610         }
1611         if (!do_ene || !do_log)
1612         {
1613             /* Write final energy file entries */
1614             print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1615                        !do_log ? fplog : nullptr, step, step, eprNORMAL,
1616                        mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1617         }
1618     }
1619
1620     /* Print some stuff... */
1621     if (MASTER(cr))
1622     {
1623         fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1624     }
1625
1626     /* IMPORTANT!
1627      * For accurate normal mode calculation it is imperative that we
1628      * store the last conformation into the full precision binary trajectory.
1629      *
1630      * However, we should only do it if we did NOT already write this step
1631      * above (which we did if do_x or do_f was true).
1632      */
1633     do_x = !do_per_step(step, inputrec->nstxout);
1634     do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1635
1636     write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1637                   top_global, inputrec, step,
1638                   s_min, state_global, observablesHistory);
1639
1640
1641     if (MASTER(cr))
1642     {
1643         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1644         print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1645                         s_min, sqrtNumAtoms);
1646         print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1647                         s_min, sqrtNumAtoms);
1648
1649         fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1650     }
1651
1652     finish_em(cr, outf, walltime_accounting, wcycle);
1653
1654     /* To print the actual number of steps we needed somewhere */
1655     walltime_accounting_set_nsteps_done(walltime_accounting, step);
1656 }
1657
1658
1659 void
1660 Integrator::do_lbfgs()
1661 {
1662     static const char *LBFGS = "Low-Memory BFGS Minimizer";
1663     em_state_t         ems;
1664     gmx_localtop_t    *top;
1665     gmx_enerdata_t    *enerd;
1666     gmx_global_stat_t  gstat;
1667     t_graph           *graph;
1668     int                ncorr, nmaxcorr, point, cp, neval, nminstep;
1669     double             stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1670     real              *rho, *alpha, *p, *s, **dx, **dg;
1671     real               a, b, c, maxdelta, delta;
1672     real               diag, Epot0;
1673     real               dgdx, dgdg, sq, yr, beta;
1674     t_mdebin          *mdebin;
1675     gmx_bool           converged;
1676     rvec               mu_tot;
1677     gmx_bool           do_log, do_ene, do_x, do_f, foundlower, *frozen;
1678     tensor             vir, pres;
1679     int                start, end, number_steps;
1680     gmx_mdoutf_t       outf;
1681     int                i, k, m, n, gf, step;
1682     int                mdof_flags;
1683     auto               mdatoms = mdAtoms->mdatoms();
1684
1685     if (PAR(cr))
1686     {
1687         gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1688     }
1689
1690     if (nullptr != constr)
1691     {
1692         gmx_fatal(FARGS, "The combination of constraints and L-BFGS minimization is not implemented. Either do not use constraints, or use another minimizer (e.g. steepest descent).");
1693     }
1694
1695     n        = 3*state_global->natoms;
1696     nmaxcorr = inputrec->nbfgscorr;
1697
1698     snew(frozen, n);
1699
1700     snew(p, n);
1701     snew(rho, nmaxcorr);
1702     snew(alpha, nmaxcorr);
1703
1704     snew(dx, nmaxcorr);
1705     for (i = 0; i < nmaxcorr; i++)
1706     {
1707         snew(dx[i], n);
1708     }
1709
1710     snew(dg, nmaxcorr);
1711     for (i = 0; i < nmaxcorr; i++)
1712     {
1713         snew(dg[i], n);
1714     }
1715
1716     step  = 0;
1717     neval = 0;
1718
1719     /* Init em */
1720     init_em(fplog, LBFGS, cr, ms, outputProvider, inputrec, mdrunOptions,
1721             state_global, top_global, &ems, &top,
1722             nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
1723             vsite, constr, nullptr,
1724             nfile, fnm, &outf, &mdebin, wcycle);
1725
1726     start = 0;
1727     end   = mdatoms->homenr;
1728
1729     /* We need 4 working states */
1730     em_state_t  s0 {}, s1 {}, s2 {}, s3 {};
1731     em_state_t *sa   = &s0;
1732     em_state_t *sb   = &s1;
1733     em_state_t *sc   = &s2;
1734     em_state_t *last = &s3;
1735     /* Initialize by copying the state from ems (we could skip x and f here) */
1736     *sa              = ems;
1737     *sb              = ems;
1738     *sc              = ems;
1739
1740     /* Print to log file */
1741     print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1742
1743     do_log = do_ene = do_x = do_f = TRUE;
1744
1745     /* Max number of steps */
1746     number_steps = inputrec->nsteps;
1747
1748     /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1749     gf = 0;
1750     for (i = start; i < end; i++)
1751     {
1752         if (mdatoms->cFREEZE)
1753         {
1754             gf = mdatoms->cFREEZE[i];
1755         }
1756         for (m = 0; m < DIM; m++)
1757         {
1758             frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1759         }
1760     }
1761     if (MASTER(cr))
1762     {
1763         sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1764     }
1765     if (fplog)
1766     {
1767         sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1768     }
1769
1770     if (vsite)
1771     {
1772         construct_vsites(vsite, as_rvec_array(state_global->x.data()), 1, nullptr,
1773                          top->idef.iparams, top->idef.il,
1774                          fr->ePBC, fr->bMolPBC, cr, state_global->box);
1775     }
1776
1777     /* Call the force routine and some auxiliary (neighboursearching etc.) */
1778     /* do_force always puts the charge groups in the box and shifts again
1779      * We do not unshift, so molecules are always whole
1780      */
1781     neval++;
1782     EnergyEvaluator energyEvaluator {
1783         fplog, cr, ms,
1784         top_global, top,
1785         inputrec, nrnb, wcycle, gstat,
1786         vsite, constr, fcd, graph,
1787         mdAtoms, fr, enerd
1788     };
1789     energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE);
1790
1791     if (MASTER(cr))
1792     {
1793         /* Copy stuff to the energy bin for easy printing etc. */
1794         upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1795                    mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
1796                    nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1797
1798         print_ebin_header(fplog, step, step);
1799         print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1800                    mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1801     }
1802
1803     /* Set the initial step.
1804      * since it will be multiplied by the non-normalized search direction
1805      * vector (force vector the first time), we scale it by the
1806      * norm of the force.
1807      */
1808
1809     if (MASTER(cr))
1810     {
1811         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1812         fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1813         fprintf(stderr, "   F-max             = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1814         fprintf(stderr, "   F-Norm            = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1815         fprintf(stderr, "\n");
1816         /* and copy to the log file too... */
1817         fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1818         fprintf(fplog, "   F-max             = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1819         fprintf(fplog, "   F-Norm            = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1820         fprintf(fplog, "\n");
1821     }
1822
1823     // Point is an index to the memory of search directions, where 0 is the first one.
1824     point = 0;
1825
1826     // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1827     real *fInit = static_cast<real *>(as_rvec_array(ems.f.data())[0]);
1828     for (i = 0; i < n; i++)
1829     {
1830         if (!frozen[i])
1831         {
1832             dx[point][i] = fInit[i]; /* Initial search direction */
1833         }
1834         else
1835         {
1836             dx[point][i] = 0;
1837         }
1838     }
1839
1840     // Stepsize will be modified during the search, and actually it is not critical
1841     // (the main efficiency in the algorithm comes from changing directions), but
1842     // we still need an initial value, so estimate it as the inverse of the norm
1843     // so we take small steps where the potential fluctuates a lot.
1844     stepsize  = 1.0/ems.fnorm;
1845
1846     /* Start the loop over BFGS steps.
1847      * Each successful step is counted, and we continue until
1848      * we either converge or reach the max number of steps.
1849      */
1850
1851     ncorr = 0;
1852
1853     /* Set the gradient from the force */
1854     converged = FALSE;
1855     for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1856     {
1857
1858         /* Write coordinates if necessary */
1859         do_x = do_per_step(step, inputrec->nstxout);
1860         do_f = do_per_step(step, inputrec->nstfout);
1861
1862         mdof_flags = 0;
1863         if (do_x)
1864         {
1865             mdof_flags |= MDOF_X;
1866         }
1867
1868         if (do_f)
1869         {
1870             mdof_flags |= MDOF_F;
1871         }
1872
1873         if (inputrec->bIMD)
1874         {
1875             mdof_flags |= MDOF_IMD;
1876         }
1877
1878         mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1879                                          top_global, step, (real)step, &ems.s, state_global, observablesHistory, ems.f);
1880
1881         /* Do the linesearching in the direction dx[point][0..(n-1)] */
1882
1883         /* make s a pointer to current search direction - point=0 first time we get here */
1884         s = dx[point];
1885
1886         real *xx = static_cast<real *>(as_rvec_array(ems.s.x.data())[0]);
1887         real *ff = static_cast<real *>(as_rvec_array(ems.f.data())[0]);
1888
1889         // calculate line gradient in position A
1890         for (gpa = 0, i = 0; i < n; i++)
1891         {
1892             gpa -= s[i]*ff[i];
1893         }
1894
1895         /* Calculate minimum allowed stepsize along the line, before the average (norm)
1896          * relative change in coordinate is smaller than precision
1897          */
1898         for (minstep = 0, i = 0; i < n; i++)
1899         {
1900             tmp = fabs(xx[i]);
1901             if (tmp < 1.0)
1902             {
1903                 tmp = 1.0;
1904             }
1905             tmp      = s[i]/tmp;
1906             minstep += tmp*tmp;
1907         }
1908         minstep = GMX_REAL_EPS/sqrt(minstep/n);
1909
1910         if (stepsize < minstep)
1911         {
1912             converged = TRUE;
1913             break;
1914         }
1915
1916         // Before taking any steps along the line, store the old position
1917         *last       = ems;
1918         real *lastx = static_cast<real *>(as_rvec_array(last->s.x.data())[0]);
1919         real *lastf = static_cast<real *>(as_rvec_array(last->f.data())[0]);
1920         Epot0       = ems.epot;
1921
1922         *sa         = ems;
1923
1924         /* Take a step downhill.
1925          * In theory, we should find the actual minimum of the function in this
1926          * direction, somewhere along the line.
1927          * That is quite possible, but it turns out to take 5-10 function evaluations
1928          * for each line. However, we dont really need to find the exact minimum -
1929          * it is much better to start a new BFGS step in a modified direction as soon
1930          * as we are close to it. This will save a lot of energy evaluations.
1931          *
1932          * In practice, we just try to take a single step.
1933          * If it worked (i.e. lowered the energy), we increase the stepsize but
1934          * continue straight to the next BFGS step without trying to find any minimum,
1935          * i.e. we change the search direction too. If the line was smooth, it is
1936          * likely we are in a smooth region, and then it makes sense to take longer
1937          * steps in the modified search direction too.
1938          *
1939          * If it didn't work (higher energy), there must be a minimum somewhere between
1940          * the old position and the new one. Then we need to start by finding a lower
1941          * value before we change search direction. Since the energy was apparently
1942          * quite rough, we need to decrease the step size.
1943          *
1944          * Due to the finite numerical accuracy, it turns out that it is a good idea
1945          * to accept a SMALL increase in energy, if the derivative is still downhill.
1946          * This leads to lower final energies in the tests I've done. / Erik
1947          */
1948
1949         // State "A" is the first position along the line.
1950         // reference position along line is initially zero
1951         a          = 0.0;
1952
1953         // Check stepsize first. We do not allow displacements
1954         // larger than emstep.
1955         //
1956         do
1957         {
1958             // Pick a new position C by adding stepsize to A.
1959             c        = a + stepsize;
1960
1961             // Calculate what the largest change in any individual coordinate
1962             // would be (translation along line * gradient along line)
1963             maxdelta = 0;
1964             for (i = 0; i < n; i++)
1965             {
1966                 delta = c*s[i];
1967                 if (delta > maxdelta)
1968                 {
1969                     maxdelta = delta;
1970                 }
1971             }
1972             // If any displacement is larger than the stepsize limit, reduce the step
1973             if (maxdelta > inputrec->em_stepsize)
1974             {
1975                 stepsize *= 0.1;
1976             }
1977         }
1978         while (maxdelta > inputrec->em_stepsize);
1979
1980         // Take a trial step and move the coordinate array xc[] to position C
1981         real *xc = static_cast<real *>(as_rvec_array(sc->s.x.data())[0]);
1982         for (i = 0; i < n; i++)
1983         {
1984             xc[i] = lastx[i] + c*s[i];
1985         }
1986
1987         neval++;
1988         // Calculate energy for the trial step in position C
1989         energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE);
1990
1991         // Calc line gradient in position C
1992         real *fc = static_cast<real *>(as_rvec_array(sc->f.data())[0]);
1993         for (gpc = 0, i = 0; i < n; i++)
1994         {
1995             gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1996         }
1997         /* Sum the gradient along the line across CPUs */
1998         if (PAR(cr))
1999         {
2000             gmx_sumd(1, &gpc, cr);
2001         }
2002
2003         // This is the max amount of increase in energy we tolerate.
2004         // By allowing VERY small changes (close to numerical precision) we
2005         // frequently find even better (lower) final energies.
2006         tmp = std::sqrt(GMX_REAL_EPS)*fabs(sa->epot);
2007
2008         // Accept the step if the energy is lower in the new position C (compared to A),
2009         // or if it is not significantly higher and the line derivative is still negative.
2010         if (sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp)))
2011         {
2012             // Great, we found a better energy. We no longer try to alter the
2013             // stepsize, but simply accept this new better position. The we select a new
2014             // search direction instead, which will be much more efficient than continuing
2015             // to take smaller steps along a line. Set fnorm based on the new C position,
2016             // which will be used to update the stepsize to 1/fnorm further down.
2017             foundlower = TRUE;
2018         }
2019         else
2020         {
2021             // If we got here, the energy is NOT lower in point C, i.e. it will be the same
2022             // or higher than in point A. In this case it is pointless to move to point C,
2023             // so we will have to do more iterations along the same line to find a smaller
2024             // value in the interval [A=0.0,C].
2025             // Here, A is still 0.0, but that will change when we do a search in the interval
2026             // [0.0,C] below. That search we will do by interpolation or bisection rather
2027             // than with the stepsize, so no need to modify it. For the next search direction
2028             // it will be reset to 1/fnorm anyway.
2029             foundlower = FALSE;
2030         }
2031
2032         if (!foundlower)
2033         {
2034             // OK, if we didn't find a lower value we will have to locate one now - there must
2035             // be one in the interval [a,c].
2036             // The same thing is valid here, though: Don't spend dozens of iterations to find
2037             // the line minimum. We try to interpolate based on the derivative at the endpoints,
2038             // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2039             // I also have a safeguard for potentially really pathological functions so we never
2040             // take more than 20 steps before we give up.
2041             // If we already found a lower value we just skip this step and continue to the update.
2042             real fnorm = 0;
2043             nminstep   = 0;
2044             do
2045             {
2046                 // Select a new trial point B in the interval [A,C].
2047                 // If the derivatives at points a & c have different sign we interpolate to zero,
2048                 // otherwise just do a bisection since there might be multiple minima/maxima
2049                 // inside the interval.
2050                 if (gpa < 0 && gpc > 0)
2051                 {
2052                     b = a + gpa*(a-c)/(gpc-gpa);
2053                 }
2054                 else
2055                 {
2056                     b = 0.5*(a+c);
2057                 }
2058
2059                 /* safeguard if interpolation close to machine accuracy causes errors:
2060                  * never go outside the interval
2061                  */
2062                 if (b <= a || b >= c)
2063                 {
2064                     b = 0.5*(a+c);
2065                 }
2066
2067                 // Take a trial step to point B
2068                 real *xb = static_cast<real *>(as_rvec_array(sb->s.x.data())[0]);
2069                 for (i = 0; i < n; i++)
2070                 {
2071                     xb[i] = lastx[i] + b*s[i];
2072                 }
2073
2074                 neval++;
2075                 // Calculate energy for the trial step in point B
2076                 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE);
2077                 fnorm = sb->fnorm;
2078
2079                 // Calculate gradient in point B
2080                 real *fb = static_cast<real *>(as_rvec_array(sb->f.data())[0]);
2081                 for (gpb = 0, i = 0; i < n; i++)
2082                 {
2083                     gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2084
2085                 }
2086                 /* Sum the gradient along the line across CPUs */
2087                 if (PAR(cr))
2088                 {
2089                     gmx_sumd(1, &gpb, cr);
2090                 }
2091
2092                 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2093                 // at the new point B, and rename the endpoints of this new interval A and C.
2094                 if (gpb > 0)
2095                 {
2096                     /* Replace c endpoint with b */
2097                     c   = b;
2098                     /* swap states b and c */
2099                     swap_em_state(&sb, &sc);
2100                 }
2101                 else
2102                 {
2103                     /* Replace a endpoint with b */
2104                     a   = b;
2105                     /* swap states a and b */
2106                     swap_em_state(&sa, &sb);
2107                 }
2108
2109                 /*
2110                  * Stop search as soon as we find a value smaller than the endpoints,
2111                  * or if the tolerance is below machine precision.
2112                  * Never run more than 20 steps, no matter what.
2113                  */
2114                 nminstep++;
2115             }
2116             while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2117
2118             if (std::fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2119             {
2120                 /* OK. We couldn't find a significantly lower energy.
2121                  * If ncorr==0 this was steepest descent, and then we give up.
2122                  * If not, reset memory to restart as steepest descent before quitting.
2123                  */
2124                 if (ncorr == 0)
2125                 {
2126                     /* Converged */
2127                     converged = TRUE;
2128                     break;
2129                 }
2130                 else
2131                 {
2132                     /* Reset memory */
2133                     ncorr = 0;
2134                     /* Search in gradient direction */
2135                     for (i = 0; i < n; i++)
2136                     {
2137                         dx[point][i] = ff[i];
2138                     }
2139                     /* Reset stepsize */
2140                     stepsize = 1.0/fnorm;
2141                     continue;
2142                 }
2143             }
2144
2145             /* Select min energy state of A & C, put the best in xx/ff/Epot
2146              */
2147             if (sc->epot < sa->epot)
2148             {
2149                 /* Use state C */
2150                 ems        = *sc;
2151                 step_taken = c;
2152             }
2153             else
2154             {
2155                 /* Use state A */
2156                 ems        = *sa;
2157                 step_taken = a;
2158             }
2159
2160         }
2161         else
2162         {
2163             /* found lower */
2164             /* Use state C */
2165             ems        = *sc;
2166             step_taken = c;
2167         }
2168
2169         /* Update the memory information, and calculate a new
2170          * approximation of the inverse hessian
2171          */
2172
2173         /* Have new data in Epot, xx, ff */
2174         if (ncorr < nmaxcorr)
2175         {
2176             ncorr++;
2177         }
2178
2179         for (i = 0; i < n; i++)
2180         {
2181             dg[point][i]  = lastf[i]-ff[i];
2182             dx[point][i] *= step_taken;
2183         }
2184
2185         dgdg = 0;
2186         dgdx = 0;
2187         for (i = 0; i < n; i++)
2188         {
2189             dgdg += dg[point][i]*dg[point][i];
2190             dgdx += dg[point][i]*dx[point][i];
2191         }
2192
2193         diag = dgdx/dgdg;
2194
2195         rho[point] = 1.0/dgdx;
2196         point++;
2197
2198         if (point >= nmaxcorr)
2199         {
2200             point = 0;
2201         }
2202
2203         /* Update */
2204         for (i = 0; i < n; i++)
2205         {
2206             p[i] = ff[i];
2207         }
2208
2209         cp = point;
2210
2211         /* Recursive update. First go back over the memory points */
2212         for (k = 0; k < ncorr; k++)
2213         {
2214             cp--;
2215             if (cp < 0)
2216             {
2217                 cp = ncorr-1;
2218             }
2219
2220             sq = 0;
2221             for (i = 0; i < n; i++)
2222             {
2223                 sq += dx[cp][i]*p[i];
2224             }
2225
2226             alpha[cp] = rho[cp]*sq;
2227
2228             for (i = 0; i < n; i++)
2229             {
2230                 p[i] -= alpha[cp]*dg[cp][i];
2231             }
2232         }
2233
2234         for (i = 0; i < n; i++)
2235         {
2236             p[i] *= diag;
2237         }
2238
2239         /* And then go forward again */
2240         for (k = 0; k < ncorr; k++)
2241         {
2242             yr = 0;
2243             for (i = 0; i < n; i++)
2244             {
2245                 yr += p[i]*dg[cp][i];
2246             }
2247
2248             beta = rho[cp]*yr;
2249             beta = alpha[cp]-beta;
2250
2251             for (i = 0; i < n; i++)
2252             {
2253                 p[i] += beta*dx[cp][i];
2254             }
2255
2256             cp++;
2257             if (cp >= ncorr)
2258             {
2259                 cp = 0;
2260             }
2261         }
2262
2263         for (i = 0; i < n; i++)
2264         {
2265             if (!frozen[i])
2266             {
2267                 dx[point][i] = p[i];
2268             }
2269             else
2270             {
2271                 dx[point][i] = 0;
2272             }
2273         }
2274
2275         /* Print it if necessary */
2276         if (MASTER(cr))
2277         {
2278             if (mdrunOptions.verbose)
2279             {
2280                 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2281                 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2282                         step, ems.epot, ems.fnorm/sqrtNumAtoms, ems.fmax, ems.a_fmax + 1);
2283                 fflush(stderr);
2284             }
2285             /* Store the new (lower) energies */
2286             upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2287                        mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
2288                        nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2289             do_log = do_per_step(step, inputrec->nstlog);
2290             do_ene = do_per_step(step, inputrec->nstenergy);
2291             if (do_log)
2292             {
2293                 print_ebin_header(fplog, step, step);
2294             }
2295             print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2296                        do_log ? fplog : nullptr, step, step, eprNORMAL,
2297                        mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
2298         }
2299
2300         /* Send x and E to IMD client, if bIMD is TRUE. */
2301         if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, as_rvec_array(state_global->x.data()), inputrec, 0, wcycle) && MASTER(cr))
2302         {
2303             IMD_send_positions(inputrec->imd);
2304         }
2305
2306         // Reset stepsize in we are doing more iterations
2307         stepsize = 1.0/ems.fnorm;
2308
2309         /* Stop when the maximum force lies below tolerance.
2310          * If we have reached machine precision, converged is already set to true.
2311          */
2312         converged = converged || (ems.fmax < inputrec->em_tol);
2313
2314     }   /* End of the loop */
2315
2316     /* IMD cleanup, if bIMD is TRUE. */
2317     IMD_finalize(inputrec->bIMD, inputrec->imd);
2318
2319     if (converged)
2320     {
2321         step--; /* we never took that last step in this case */
2322
2323     }
2324     if (ems.fmax > inputrec->em_tol)
2325     {
2326         if (MASTER(cr))
2327         {
2328             warn_step(fplog, inputrec->em_tol, ems.fmax,
2329                       step-1 == number_steps, FALSE);
2330         }
2331         converged = FALSE;
2332     }
2333
2334     /* If we printed energy and/or logfile last step (which was the last step)
2335      * we don't have to do it again, but otherwise print the final values.
2336      */
2337     if (!do_log) /* Write final value to log since we didn't do anythin last step */
2338     {
2339         print_ebin_header(fplog, step, step);
2340     }
2341     if (!do_ene || !do_log) /* Write final energy file entries */
2342     {
2343         print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2344                    !do_log ? fplog : nullptr, step, step, eprNORMAL,
2345                    mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
2346     }
2347
2348     /* Print some stuff... */
2349     if (MASTER(cr))
2350     {
2351         fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2352     }
2353
2354     /* IMPORTANT!
2355      * For accurate normal mode calculation it is imperative that we
2356      * store the last conformation into the full precision binary trajectory.
2357      *
2358      * However, we should only do it if we did NOT already write this step
2359      * above (which we did if do_x or do_f was true).
2360      */
2361     do_x = !do_per_step(step, inputrec->nstxout);
2362     do_f = !do_per_step(step, inputrec->nstfout);
2363     write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2364                   top_global, inputrec, step,
2365                   &ems, state_global, observablesHistory);
2366
2367     if (MASTER(cr))
2368     {
2369         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2370         print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2371                         number_steps, &ems, sqrtNumAtoms);
2372         print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2373                         number_steps, &ems, sqrtNumAtoms);
2374
2375         fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2376     }
2377
2378     finish_em(cr, outf, walltime_accounting, wcycle);
2379
2380     /* To print the actual number of steps we needed somewhere */
2381     walltime_accounting_set_nsteps_done(walltime_accounting, step);
2382 }
2383
2384 void
2385 Integrator::do_steep()
2386 {
2387     const char       *SD = "Steepest Descents";
2388     gmx_localtop_t   *top;
2389     gmx_enerdata_t   *enerd;
2390     gmx_global_stat_t gstat;
2391     t_graph          *graph;
2392     real              stepsize;
2393     real              ustep;
2394     gmx_mdoutf_t      outf;
2395     t_mdebin         *mdebin;
2396     gmx_bool          bDone, bAbort, do_x, do_f;
2397     tensor            vir, pres;
2398     rvec              mu_tot;
2399     int               nsteps;
2400     int               count          = 0;
2401     int               steps_accepted = 0;
2402     auto              mdatoms        = mdAtoms->mdatoms();
2403
2404     /* Create 2 states on the stack and extract pointers that we will swap */
2405     em_state_t  s0 {}, s1 {};
2406     em_state_t *s_min = &s0;
2407     em_state_t *s_try = &s1;
2408
2409     /* Init em and store the local state in s_try */
2410     init_em(fplog, SD, cr, ms, outputProvider, inputrec, mdrunOptions,
2411             state_global, top_global, s_try, &top,
2412             nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
2413             vsite, constr, nullptr,
2414             nfile, fnm, &outf, &mdebin, wcycle);
2415
2416     /* Print to log file  */
2417     print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2418
2419     /* Set variables for stepsize (in nm). This is the largest
2420      * step that we are going to make in any direction.
2421      */
2422     ustep    = inputrec->em_stepsize;
2423     stepsize = 0;
2424
2425     /* Max number of steps  */
2426     nsteps = inputrec->nsteps;
2427
2428     if (MASTER(cr))
2429     {
2430         /* Print to the screen  */
2431         sp_header(stderr, SD, inputrec->em_tol, nsteps);
2432     }
2433     if (fplog)
2434     {
2435         sp_header(fplog, SD, inputrec->em_tol, nsteps);
2436     }
2437     EnergyEvaluator energyEvaluator {
2438         fplog, cr, ms,
2439         top_global, top,
2440         inputrec, nrnb, wcycle, gstat,
2441         vsite, constr, fcd, graph,
2442         mdAtoms, fr, enerd
2443     };
2444
2445     /**** HERE STARTS THE LOOP ****
2446      * count is the counter for the number of steps
2447      * bDone will be TRUE when the minimization has converged
2448      * bAbort will be TRUE when nsteps steps have been performed or when
2449      * the stepsize becomes smaller than is reasonable for machine precision
2450      */
2451     count  = 0;
2452     bDone  = FALSE;
2453     bAbort = FALSE;
2454     while (!bDone && !bAbort)
2455     {
2456         bAbort = (nsteps >= 0) && (count == nsteps);
2457
2458         /* set new coordinates, except for first step */
2459         bool validStep = true;
2460         if (count > 0)
2461         {
2462             validStep =
2463                 do_em_step(cr, inputrec, mdatoms,
2464                            s_min, stepsize, &s_min->f, s_try,
2465                            constr, count);
2466         }
2467
2468         if (validStep)
2469         {
2470             energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0);
2471         }
2472         else
2473         {
2474             // Signal constraint error during stepping with energy=inf
2475             s_try->epot = std::numeric_limits<real>::infinity();
2476         }
2477
2478         if (MASTER(cr))
2479         {
2480             print_ebin_header(fplog, count, count);
2481         }
2482
2483         if (count == 0)
2484         {
2485             s_min->epot = s_try->epot;
2486         }
2487
2488         /* Print it if necessary  */
2489         if (MASTER(cr))
2490         {
2491             if (mdrunOptions.verbose)
2492             {
2493                 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2494                         count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2495                         ( (count == 0) || (s_try->epot < s_min->epot) ) ? '\n' : '\r');
2496                 fflush(stderr);
2497             }
2498
2499             if ( (count == 0) || (s_try->epot < s_min->epot) )
2500             {
2501                 /* Store the new (lower) energies  */
2502                 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2503                            mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2504                            s_try->s.box, nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2505
2506                 /* Prepare IMD energy record, if bIMD is TRUE. */
2507                 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2508
2509                 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2510                            do_per_step(steps_accepted, inputrec->nstdisreout),
2511                            do_per_step(steps_accepted, inputrec->nstorireout),
2512                            fplog, count, count, eprNORMAL,
2513                            mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
2514                 fflush(fplog);
2515             }
2516         }
2517
2518         /* Now if the new energy is smaller than the previous...
2519          * or if this is the first step!
2520          * or if we did random steps!
2521          */
2522
2523         if ( (count == 0) || (s_try->epot < s_min->epot) )
2524         {
2525             steps_accepted++;
2526
2527             /* Test whether the convergence criterion is met...  */
2528             bDone = (s_try->fmax < inputrec->em_tol);
2529
2530             /* Copy the arrays for force, positions and energy  */
2531             /* The 'Min' array always holds the coords and forces of the minimal
2532                sampled energy  */
2533             swap_em_state(&s_min, &s_try);
2534             if (count > 0)
2535             {
2536                 ustep *= 1.2;
2537             }
2538
2539             /* Write to trn, if necessary */
2540             do_x = do_per_step(steps_accepted, inputrec->nstxout);
2541             do_f = do_per_step(steps_accepted, inputrec->nstfout);
2542             write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
2543                           top_global, inputrec, count,
2544                           s_min, state_global, observablesHistory);
2545         }
2546         else
2547         {
2548             /* If energy is not smaller make the step smaller...  */
2549             ustep *= 0.5;
2550
2551             if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2552             {
2553                 /* Reload the old state */
2554                 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2555                                        s_min, top, mdAtoms, fr, vsite, constr,
2556                                        nrnb, wcycle);
2557             }
2558         }
2559
2560         /* Determine new step  */
2561         stepsize = ustep/s_min->fmax;
2562
2563         /* Check if stepsize is too small, with 1 nm as a characteristic length */
2564 #if GMX_DOUBLE
2565         if (count == nsteps || ustep < 1e-12)
2566 #else
2567         if (count == nsteps || ustep < 1e-6)
2568 #endif
2569         {
2570             if (MASTER(cr))
2571             {
2572                 warn_step(fplog, inputrec->em_tol, s_min->fmax,
2573                           count == nsteps, constr != nullptr);
2574             }
2575             bAbort = TRUE;
2576         }
2577
2578         /* Send IMD energies and positions, if bIMD is TRUE. */
2579         if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box,
2580                    MASTER(cr) ? as_rvec_array(state_global->x.data()) : nullptr,
2581                    inputrec, 0, wcycle) &&
2582             MASTER(cr))
2583         {
2584             IMD_send_positions(inputrec->imd);
2585         }
2586
2587         count++;
2588     }   /* End of the loop  */
2589
2590     /* IMD cleanup, if bIMD is TRUE. */
2591     IMD_finalize(inputrec->bIMD, inputrec->imd);
2592
2593     /* Print some data...  */
2594     if (MASTER(cr))
2595     {
2596         fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2597     }
2598     write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2599                   top_global, inputrec, count,
2600                   s_min, state_global, observablesHistory);
2601
2602     if (MASTER(cr))
2603     {
2604         double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2605
2606         print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2607                         s_min, sqrtNumAtoms);
2608         print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2609                         s_min, sqrtNumAtoms);
2610     }
2611
2612     finish_em(cr, outf, walltime_accounting, wcycle);
2613
2614     /* To print the actual number of steps we needed somewhere */
2615     inputrec->nsteps = count;
2616
2617     walltime_accounting_set_nsteps_done(walltime_accounting, count);
2618 }
2619
2620 void
2621 Integrator::do_nm()
2622 {
2623     const char          *NM = "Normal Mode Analysis";
2624     gmx_mdoutf_t         outf;
2625     int                  nnodes, node;
2626     gmx_localtop_t      *top;
2627     gmx_enerdata_t      *enerd;
2628     gmx_global_stat_t    gstat;
2629     t_graph             *graph;
2630     tensor               vir, pres;
2631     rvec                 mu_tot;
2632     rvec                *fneg, *dfdx;
2633     gmx_bool             bSparse; /* use sparse matrix storage format */
2634     size_t               sz;
2635     gmx_sparsematrix_t * sparse_matrix           = nullptr;
2636     real           *     full_matrix             = nullptr;
2637
2638     /* added with respect to mdrun */
2639     int                       row, col;
2640     real                      der_range = 10.0*std::sqrt(GMX_REAL_EPS);
2641     real                      x_min;
2642     bool                      bIsMaster = MASTER(cr);
2643     auto                      mdatoms   = mdAtoms->mdatoms();
2644
2645     if (constr != nullptr)
2646     {
2647         gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2648     }
2649
2650     gmx_shellfc_t *shellfc;
2651
2652     em_state_t     state_work {};
2653
2654     /* Init em and store the local state in state_minimum */
2655     init_em(fplog, NM, cr, ms, outputProvider, inputrec, mdrunOptions,
2656             state_global, top_global, &state_work, &top,
2657             nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
2658             vsite, constr, &shellfc,
2659             nfile, fnm, &outf, nullptr, wcycle);
2660
2661     std::vector<size_t> atom_index = get_atom_index(top_global);
2662     snew(fneg, atom_index.size());
2663     snew(dfdx, atom_index.size());
2664
2665 #if !GMX_DOUBLE
2666     if (bIsMaster)
2667     {
2668         fprintf(stderr,
2669                 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2670                 "      which MIGHT not be accurate enough for normal mode analysis.\n"
2671                 "      GROMACS now uses sparse matrix storage, so the memory requirements\n"
2672                 "      are fairly modest even if you recompile in double precision.\n\n");
2673     }
2674 #endif
2675
2676     /* Check if we can/should use sparse storage format.
2677      *
2678      * Sparse format is only useful when the Hessian itself is sparse, which it
2679      * will be when we use a cutoff.
2680      * For small systems (n<1000) it is easier to always use full matrix format, though.
2681      */
2682     if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
2683     {
2684         GMX_LOG(mdlog.warning).appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2685         bSparse = FALSE;
2686     }
2687     else if (atom_index.size() < 1000)
2688     {
2689         GMX_LOG(mdlog.warning).appendTextFormatted("Small system size (N=%d), using full Hessian format.",
2690                                                    atom_index.size());
2691         bSparse = FALSE;
2692     }
2693     else
2694     {
2695         GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2696         bSparse = TRUE;
2697     }
2698
2699     /* Number of dimensions, based on real atoms, that is not vsites or shell */
2700     sz = DIM*atom_index.size();
2701
2702     fprintf(stderr, "Allocating Hessian memory...\n\n");
2703
2704     if (bSparse)
2705     {
2706         sparse_matrix = gmx_sparsematrix_init(sz);
2707         sparse_matrix->compressed_symmetric = TRUE;
2708     }
2709     else
2710     {
2711         snew(full_matrix, sz*sz);
2712     }
2713
2714     init_nrnb(nrnb);
2715
2716
2717     /* Write start time and temperature */
2718     print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2719
2720     /* fudge nr of steps to nr of atoms */
2721     inputrec->nsteps = atom_index.size()*2;
2722
2723     if (bIsMaster)
2724     {
2725         fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2726                 *(top_global->name), (int)inputrec->nsteps);
2727     }
2728
2729     nnodes = cr->nnodes;
2730
2731     /* Make evaluate_energy do a single node force calculation */
2732     cr->nnodes = 1;
2733     EnergyEvaluator energyEvaluator {
2734         fplog, cr, ms,
2735         top_global, top,
2736         inputrec, nrnb, wcycle, gstat,
2737         vsite, constr, fcd, graph,
2738         mdAtoms, fr, enerd
2739     };
2740     energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE);
2741     cr->nnodes = nnodes;
2742
2743     /* if forces are not small, warn user */
2744     get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
2745
2746     GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
2747     if (state_work.fmax > 1.0e-3)
2748     {
2749         GMX_LOG(mdlog.warning).appendText(
2750                 "The force is probably not small enough to "
2751                 "ensure that you are at a minimum.\n"
2752                 "Be aware that negative eigenvalues may occur\n"
2753                 "when the resulting matrix is diagonalized.");
2754     }
2755
2756     /***********************************************************
2757      *
2758      *      Loop over all pairs in matrix
2759      *
2760      *      do_force called twice. Once with positive and
2761      *      once with negative displacement
2762      *
2763      ************************************************************/
2764
2765     /* Steps are divided one by one over the nodes */
2766     bool bNS = true;
2767     for (unsigned int aid = cr->nodeid; aid < atom_index.size(); aid += nnodes)
2768     {
2769         size_t atom = atom_index[aid];
2770         for (size_t d = 0; d < DIM; d++)
2771         {
2772             gmx_int64_t step        = 0;
2773             int         force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2774             double      t           = 0;
2775
2776             x_min = state_work.s.x[atom][d];
2777
2778             for (unsigned int dx = 0; (dx < 2); dx++)
2779             {
2780                 if (dx == 0)
2781                 {
2782                     state_work.s.x[atom][d] = x_min - der_range;
2783                 }
2784                 else
2785                 {
2786                     state_work.s.x[atom][d] = x_min + der_range;
2787                 }
2788
2789                 /* Make evaluate_energy do a single node force calculation */
2790                 cr->nnodes = 1;
2791                 if (shellfc)
2792                 {
2793                     /* Now is the time to relax the shells */
2794                     relax_shell_flexcon(fplog,
2795                                         cr,
2796                                         ms,
2797                                         mdrunOptions.verbose,
2798                                         step,
2799                                         inputrec,
2800                                         bNS,
2801                                         force_flags,
2802                                         top,
2803                                         constr,
2804                                         enerd,
2805                                         fcd,
2806                                         &state_work.s,
2807                                         state_work.f,
2808                                         vir,
2809                                         mdatoms,
2810                                         nrnb,
2811                                         wcycle,
2812                                         graph,
2813                                         &top_global->groups,
2814                                         shellfc,
2815                                         fr,
2816                                         t,
2817                                         mu_tot,
2818                                         vsite,
2819                                         DdOpenBalanceRegionBeforeForceComputation::no,
2820                                         DdCloseBalanceRegionAfterForceComputation::no);
2821                     bNS = false;
2822                     step++;
2823                 }
2824                 else
2825                 {
2826                     energyEvaluator.run(&state_work, mu_tot, vir, pres, atom*2+dx, FALSE);
2827                 }
2828
2829                 cr->nnodes = nnodes;
2830
2831                 if (dx == 0)
2832                 {
2833                     for (size_t i = 0; i < atom_index.size(); i++)
2834                     {
2835                         copy_rvec(state_work.f[atom_index[i]], fneg[i]);
2836                     }
2837                 }
2838             }
2839
2840             /* x is restored to original */
2841             state_work.s.x[atom][d] = x_min;
2842
2843             for (size_t j = 0; j < atom_index.size(); j++)
2844             {
2845                 for (size_t k = 0; (k < DIM); k++)
2846                 {
2847                     dfdx[j][k] =
2848                         -(state_work.f[atom_index[j]][k] - fneg[j][k])/(2*der_range);
2849                 }
2850             }
2851
2852             if (!bIsMaster)
2853             {
2854 #if GMX_MPI
2855 #define mpi_type GMX_MPI_REAL
2856                 MPI_Send(dfdx[0], atom_index.size()*DIM, mpi_type, MASTER(cr),
2857                          cr->nodeid, cr->mpi_comm_mygroup);
2858 #endif
2859             }
2860             else
2861             {
2862                 for (node = 0; (node < nnodes && atom+node < atom_index.size()); node++)
2863                 {
2864                     if (node > 0)
2865                     {
2866 #if GMX_MPI
2867                         MPI_Status stat;
2868                         MPI_Recv(dfdx[0], atom_index.size()*DIM, mpi_type, node, node,
2869                                  cr->mpi_comm_mygroup, &stat);
2870 #undef mpi_type
2871 #endif
2872                     }
2873
2874                     row = (atom + node)*DIM + d;
2875
2876                     for (size_t j = 0; j < atom_index.size(); j++)
2877                     {
2878                         for (size_t k = 0; k < DIM; k++)
2879                         {
2880                             col = j*DIM + k;
2881
2882                             if (bSparse)
2883                             {
2884                                 if (col >= row && dfdx[j][k] != 0.0)
2885                                 {
2886                                     gmx_sparsematrix_increment_value(sparse_matrix,
2887                                                                      row, col, dfdx[j][k]);
2888                                 }
2889                             }
2890                             else
2891                             {
2892                                 full_matrix[row*sz+col] = dfdx[j][k];
2893                             }
2894                         }
2895                     }
2896                 }
2897             }
2898
2899             if (mdrunOptions.verbose && fplog)
2900             {
2901                 fflush(fplog);
2902             }
2903         }
2904         /* write progress */
2905         if (bIsMaster && mdrunOptions.verbose)
2906         {
2907             fprintf(stderr, "\rFinished step %d out of %d",
2908                     static_cast<int>(std::min(atom+nnodes, atom_index.size())),
2909                     static_cast<int>(atom_index.size()));
2910             fflush(stderr);
2911         }
2912     }
2913
2914     if (bIsMaster)
2915     {
2916         fprintf(stderr, "\n\nWriting Hessian...\n");
2917         gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2918     }
2919
2920     finish_em(cr, outf, walltime_accounting, wcycle);
2921
2922     walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size()*2);
2923 }
2924
2925 } // namespace gmx