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