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