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