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