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