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