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