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