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