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