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