Remove mdrun -seppot
[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(inputrec, fr, 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         enerd->term[F_DVDL_CONSTR] += dvdl_constr;
793         m_add(force_vir, shake_vir, vir);
794         wallcycle_stop(wcycle, ewcCONSTR);
795     }
796     else
797     {
798         copy_mat(force_vir, vir);
799     }
800
801     clear_mat(ekin);
802     enerd->term[F_PRES] =
803         calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
804
805     sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
806
807     if (EI_ENERGY_MINIMIZATION(inputrec->eI))
808     {
809         get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
810     }
811 }
812
813 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
814                               gmx_mtop_t *mtop,
815                               em_state_t *s_min, em_state_t *s_b)
816 {
817     rvec          *fm, *fb, *fmg;
818     t_block       *cgs_gl;
819     int            ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
820     double         partsum;
821     unsigned char *grpnrFREEZE;
822
823     if (debug)
824     {
825         fprintf(debug, "Doing reorder_partsum\n");
826     }
827
828     fm = s_min->f;
829     fb = s_b->f;
830
831     cgs_gl = dd_charge_groups_global(cr->dd);
832     index  = cgs_gl->index;
833
834     /* Collect fm in a global vector fmg.
835      * This conflicts with the spirit of domain decomposition,
836      * but to fully optimize this a much more complicated algorithm is required.
837      */
838     snew(fmg, mtop->natoms);
839
840     ncg   = s_min->s.ncg_gl;
841     cg_gl = s_min->s.cg_gl;
842     i     = 0;
843     for (c = 0; c < ncg; c++)
844     {
845         cg = cg_gl[c];
846         a0 = index[cg];
847         a1 = index[cg+1];
848         for (a = a0; a < a1; a++)
849         {
850             copy_rvec(fm[i], fmg[a]);
851             i++;
852         }
853     }
854     gmx_sum(mtop->natoms*3, fmg[0], cr);
855
856     /* Now we will determine the part of the sum for the cgs in state s_b */
857     ncg         = s_b->s.ncg_gl;
858     cg_gl       = s_b->s.cg_gl;
859     partsum     = 0;
860     i           = 0;
861     gf          = 0;
862     grpnrFREEZE = mtop->groups.grpnr[egcFREEZE];
863     for (c = 0; c < ncg; c++)
864     {
865         cg = cg_gl[c];
866         a0 = index[cg];
867         a1 = index[cg+1];
868         for (a = a0; a < a1; a++)
869         {
870             if (mdatoms->cFREEZE && grpnrFREEZE)
871             {
872                 gf = grpnrFREEZE[i];
873             }
874             for (m = 0; m < DIM; m++)
875             {
876                 if (!opts->nFreeze[gf][m])
877                 {
878                     partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
879                 }
880             }
881             i++;
882         }
883     }
884
885     sfree(fmg);
886
887     return partsum;
888 }
889
890 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
891                     gmx_mtop_t *mtop,
892                     em_state_t *s_min, em_state_t *s_b)
893 {
894     rvec  *fm, *fb;
895     double sum;
896     int    gf, i, m;
897
898     /* This is just the classical Polak-Ribiere calculation of beta;
899      * it looks a bit complicated since we take freeze groups into account,
900      * and might have to sum it in parallel runs.
901      */
902
903     if (!DOMAINDECOMP(cr) ||
904         (s_min->s.ddp_count == cr->dd->ddp_count &&
905          s_b->s.ddp_count   == cr->dd->ddp_count))
906     {
907         fm  = s_min->f;
908         fb  = s_b->f;
909         sum = 0;
910         gf  = 0;
911         /* This part of code can be incorrect with DD,
912          * since the atom ordering in s_b and s_min might differ.
913          */
914         for (i = 0; i < mdatoms->homenr; i++)
915         {
916             if (mdatoms->cFREEZE)
917             {
918                 gf = mdatoms->cFREEZE[i];
919             }
920             for (m = 0; m < DIM; m++)
921             {
922                 if (!opts->nFreeze[gf][m])
923                 {
924                     sum += (fb[i][m] - fm[i][m])*fb[i][m];
925                 }
926             }
927         }
928     }
929     else
930     {
931         /* We need to reorder cgs while summing */
932         sum = reorder_partsum(cr, opts, mdatoms, mtop, s_min, s_b);
933     }
934     if (PAR(cr))
935     {
936         gmx_sumd(1, &sum, cr);
937     }
938
939     return sum/sqr(s_min->fnorm);
940 }
941
942 double do_cg(FILE *fplog, t_commrec *cr,
943              int nfile, const t_filenm fnm[],
944              const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
945              int gmx_unused nstglobalcomm,
946              gmx_vsite_t *vsite, gmx_constr_t constr,
947              int gmx_unused stepout,
948              t_inputrec *inputrec,
949              gmx_mtop_t *top_global, t_fcdata *fcd,
950              t_state *state_global,
951              t_mdatoms *mdatoms,
952              t_nrnb *nrnb, gmx_wallcycle_t wcycle,
953              gmx_edsam_t gmx_unused ed,
954              t_forcerec *fr,
955              int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
956              gmx_membed_t gmx_unused membed,
957              real gmx_unused cpt_period, real gmx_unused max_hours,
958              const char gmx_unused *deviceOptions,
959              int imdport,
960              unsigned long gmx_unused Flags,
961              gmx_walltime_accounting_t walltime_accounting)
962 {
963     const char       *CG = "Polak-Ribiere Conjugate Gradients";
964
965     em_state_t       *s_min, *s_a, *s_b, *s_c;
966     gmx_localtop_t   *top;
967     gmx_enerdata_t   *enerd;
968     rvec             *f;
969     gmx_global_stat_t gstat;
970     t_graph          *graph;
971     rvec             *f_global, *p, *sf, *sfm;
972     double            gpa, gpb, gpc, tmp, sum[2], minstep;
973     real              fnormn;
974     real              stepsize;
975     real              a, b, c, beta = 0.0;
976     real              epot_repl = 0;
977     real              pnorm;
978     t_mdebin         *mdebin;
979     gmx_bool          converged, foundlower;
980     rvec              mu_tot;
981     gmx_bool          do_log = FALSE, do_ene = FALSE, do_x, do_f;
982     tensor            vir, pres;
983     int               number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
984     gmx_mdoutf_t      outf;
985     int               i, m, gf, step, nminstep;
986     real              terminate = 0;
987
988     step = 0;
989
990     s_min = init_em_state();
991     s_a   = init_em_state();
992     s_b   = init_em_state();
993     s_c   = init_em_state();
994
995     /* Init em and store the local state in s_min */
996     init_em(fplog, CG, cr, inputrec,
997             state_global, top_global, s_min, &top, &f, &f_global,
998             nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
999             nfile, fnm, &outf, &mdebin, imdport, Flags);
1000
1001     /* Print to log file */
1002     print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1003
1004     /* Max number of steps */
1005     number_steps = inputrec->nsteps;
1006
1007     if (MASTER(cr))
1008     {
1009         sp_header(stderr, CG, inputrec->em_tol, number_steps);
1010     }
1011     if (fplog)
1012     {
1013         sp_header(fplog, CG, inputrec->em_tol, number_steps);
1014     }
1015
1016     /* Call the force routine and some auxiliary (neighboursearching etc.) */
1017     /* do_force always puts the charge groups in the box and shifts again
1018      * We do not unshift, so molecules are always whole in congrad.c
1019      */
1020     evaluate_energy(fplog, cr,
1021                     top_global, s_min, top,
1022                     inputrec, nrnb, wcycle, gstat,
1023                     vsite, constr, fcd, graph, mdatoms, fr,
1024                     mu_tot, enerd, vir, pres, -1, TRUE);
1025     where();
1026
1027     if (MASTER(cr))
1028     {
1029         /* Copy stuff to the energy bin for easy printing etc. */
1030         upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1031                    mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1032                    NULL, NULL, vir, pres, NULL, mu_tot, constr);
1033
1034         print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1035         print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1036                    TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1037     }
1038     where();
1039
1040     /* Estimate/guess the initial stepsize */
1041     stepsize = inputrec->em_stepsize/s_min->fnorm;
1042
1043     if (MASTER(cr))
1044     {
1045         fprintf(stderr, "   F-max             = %12.5e on atom %d\n",
1046                 s_min->fmax, s_min->a_fmax+1);
1047         fprintf(stderr, "   F-Norm            = %12.5e\n",
1048                 s_min->fnorm/sqrt(state_global->natoms));
1049         fprintf(stderr, "\n");
1050         /* and copy to the log file too... */
1051         fprintf(fplog, "   F-max             = %12.5e on atom %d\n",
1052                 s_min->fmax, s_min->a_fmax+1);
1053         fprintf(fplog, "   F-Norm            = %12.5e\n",
1054                 s_min->fnorm/sqrt(state_global->natoms));
1055         fprintf(fplog, "\n");
1056     }
1057     /* Start the loop over CG steps.
1058      * Each successful step is counted, and we continue until
1059      * we either converge or reach the max number of steps.
1060      */
1061     converged = FALSE;
1062     for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1063     {
1064
1065         /* start taking steps in a new direction
1066          * First time we enter the routine, beta=0, and the direction is
1067          * simply the negative gradient.
1068          */
1069
1070         /* Calculate the new direction in p, and the gradient in this direction, gpa */
1071         p   = s_min->s.cg_p;
1072         sf  = s_min->f;
1073         gpa = 0;
1074         gf  = 0;
1075         for (i = 0; i < mdatoms->homenr; i++)
1076         {
1077             if (mdatoms->cFREEZE)
1078             {
1079                 gf = mdatoms->cFREEZE[i];
1080             }
1081             for (m = 0; m < DIM; m++)
1082             {
1083                 if (!inputrec->opts.nFreeze[gf][m])
1084                 {
1085                     p[i][m] = sf[i][m] + beta*p[i][m];
1086                     gpa    -= p[i][m]*sf[i][m];
1087                     /* f is negative gradient, thus the sign */
1088                 }
1089                 else
1090                 {
1091                     p[i][m] = 0;
1092                 }
1093             }
1094         }
1095
1096         /* Sum the gradient along the line across CPUs */
1097         if (PAR(cr))
1098         {
1099             gmx_sumd(1, &gpa, cr);
1100         }
1101
1102         /* Calculate the norm of the search vector */
1103         get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1104
1105         /* Just in case stepsize reaches zero due to numerical precision... */
1106         if (stepsize <= 0)
1107         {
1108             stepsize = inputrec->em_stepsize/pnorm;
1109         }
1110
1111         /*
1112          * Double check the value of the derivative in the search direction.
1113          * If it is positive it must be due to the old information in the
1114          * CG formula, so just remove that and start over with beta=0.
1115          * This corresponds to a steepest descent step.
1116          */
1117         if (gpa > 0)
1118         {
1119             beta = 0;
1120             step--;   /* Don't count this step since we are restarting */
1121             continue; /* Go back to the beginning of the big for-loop */
1122         }
1123
1124         /* Calculate minimum allowed stepsize, before the average (norm)
1125          * relative change in coordinate is smaller than precision
1126          */
1127         minstep = 0;
1128         for (i = 0; i < mdatoms->homenr; i++)
1129         {
1130             for (m = 0; m < DIM; m++)
1131             {
1132                 tmp = fabs(s_min->s.x[i][m]);
1133                 if (tmp < 1.0)
1134                 {
1135                     tmp = 1.0;
1136                 }
1137                 tmp      = p[i][m]/tmp;
1138                 minstep += tmp*tmp;
1139             }
1140         }
1141         /* Add up from all CPUs */
1142         if (PAR(cr))
1143         {
1144             gmx_sumd(1, &minstep, cr);
1145         }
1146
1147         minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1148
1149         if (stepsize < minstep)
1150         {
1151             converged = TRUE;
1152             break;
1153         }
1154
1155         /* Write coordinates if necessary */
1156         do_x = do_per_step(step, inputrec->nstxout);
1157         do_f = do_per_step(step, inputrec->nstfout);
1158
1159         write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1160                       top_global, inputrec, step,
1161                       s_min, state_global, f_global);
1162
1163         /* Take a step downhill.
1164          * In theory, we should minimize the function along this direction.
1165          * That is quite possible, but it turns out to take 5-10 function evaluations
1166          * for each line. However, we dont really need to find the exact minimum -
1167          * it is much better to start a new CG step in a modified direction as soon
1168          * as we are close to it. This will save a lot of energy evaluations.
1169          *
1170          * In practice, we just try to take a single step.
1171          * If it worked (i.e. lowered the energy), we increase the stepsize but
1172          * the continue straight to the next CG step without trying to find any minimum.
1173          * If it didn't work (higher energy), there must be a minimum somewhere between
1174          * the old position and the new one.
1175          *
1176          * Due to the finite numerical accuracy, it turns out that it is a good idea
1177          * to even accept a SMALL increase in energy, if the derivative is still downhill.
1178          * This leads to lower final energies in the tests I've done. / Erik
1179          */
1180         s_a->epot = s_min->epot;
1181         a         = 0.0;
1182         c         = a + stepsize; /* reference position along line is zero */
1183
1184         if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1185         {
1186             em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1187                                    s_min, top, mdatoms, fr, vsite, constr,
1188                                    nrnb, wcycle);
1189         }
1190
1191         /* Take a trial step (new coords in s_c) */
1192         do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1193                    constr, top, nrnb, wcycle, -1);
1194
1195         neval++;
1196         /* Calculate energy for the trial step */
1197         evaluate_energy(fplog, cr,
1198                         top_global, s_c, top,
1199                         inputrec, nrnb, wcycle, gstat,
1200                         vsite, constr, fcd, graph, mdatoms, fr,
1201                         mu_tot, enerd, vir, pres, -1, FALSE);
1202
1203         /* Calc derivative along line */
1204         p   = s_c->s.cg_p;
1205         sf  = s_c->f;
1206         gpc = 0;
1207         for (i = 0; i < mdatoms->homenr; i++)
1208         {
1209             for (m = 0; m < DIM; m++)
1210             {
1211                 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1212             }
1213         }
1214         /* Sum the gradient along the line across CPUs */
1215         if (PAR(cr))
1216         {
1217             gmx_sumd(1, &gpc, cr);
1218         }
1219
1220         /* This is the max amount of increase in energy we tolerate */
1221         tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1222
1223         /* Accept the step if the energy is lower, or if it is not significantly higher
1224          * and the line derivative is still negative.
1225          */
1226         if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1227         {
1228             foundlower = TRUE;
1229             /* Great, we found a better energy. Increase step for next iteration
1230              * if we are still going down, decrease it otherwise
1231              */
1232             if (gpc < 0)
1233             {
1234                 stepsize *= 1.618034; /* The golden section */
1235             }
1236             else
1237             {
1238                 stepsize *= 0.618034; /* 1/golden section */
1239             }
1240         }
1241         else
1242         {
1243             /* New energy is the same or higher. We will have to do some work
1244              * to find a smaller value in the interval. Take smaller step next time!
1245              */
1246             foundlower = FALSE;
1247             stepsize  *= 0.618034;
1248         }
1249
1250
1251
1252
1253         /* OK, if we didn't find a lower value we will have to locate one now - there must
1254          * be one in the interval [a=0,c].
1255          * The same thing is valid here, though: Don't spend dozens of iterations to find
1256          * the line minimum. We try to interpolate based on the derivative at the endpoints,
1257          * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1258          *
1259          * I also have a safeguard for potentially really patological functions so we never
1260          * take more than 20 steps before we give up ...
1261          *
1262          * If we already found a lower value we just skip this step and continue to the update.
1263          */
1264         if (!foundlower)
1265         {
1266             nminstep = 0;
1267
1268             do
1269             {
1270                 /* Select a new trial point.
1271                  * If the derivatives at points a & c have different sign we interpolate to zero,
1272                  * otherwise just do a bisection.
1273                  */
1274                 if (gpa < 0 && gpc > 0)
1275                 {
1276                     b = a + gpa*(a-c)/(gpc-gpa);
1277                 }
1278                 else
1279                 {
1280                     b = 0.5*(a+c);
1281                 }
1282
1283                 /* safeguard if interpolation close to machine accuracy causes errors:
1284                  * never go outside the interval
1285                  */
1286                 if (b <= a || b >= c)
1287                 {
1288                     b = 0.5*(a+c);
1289                 }
1290
1291                 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1292                 {
1293                     /* Reload the old state */
1294                     em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1295                                            s_min, top, mdatoms, fr, vsite, constr,
1296                                            nrnb, wcycle);
1297                 }
1298
1299                 /* Take a trial step to this new point - new coords in s_b */
1300                 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1301                            constr, top, nrnb, wcycle, -1);
1302
1303                 neval++;
1304                 /* Calculate energy for the trial step */
1305                 evaluate_energy(fplog, cr,
1306                                 top_global, s_b, top,
1307                                 inputrec, nrnb, wcycle, gstat,
1308                                 vsite, constr, fcd, graph, mdatoms, fr,
1309                                 mu_tot, enerd, vir, pres, -1, FALSE);
1310
1311                 /* p does not change within a step, but since the domain decomposition
1312                  * might change, we have to use cg_p of s_b here.
1313                  */
1314                 p   = s_b->s.cg_p;
1315                 sf  = s_b->f;
1316                 gpb = 0;
1317                 for (i = 0; i < mdatoms->homenr; i++)
1318                 {
1319                     for (m = 0; m < DIM; m++)
1320                     {
1321                         gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1322                     }
1323                 }
1324                 /* Sum the gradient along the line across CPUs */
1325                 if (PAR(cr))
1326                 {
1327                     gmx_sumd(1, &gpb, cr);
1328                 }
1329
1330                 if (debug)
1331                 {
1332                     fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1333                             s_a->epot, s_b->epot, s_c->epot, gpb);
1334                 }
1335
1336                 epot_repl = s_b->epot;
1337
1338                 /* Keep one of the intervals based on the value of the derivative at the new point */
1339                 if (gpb > 0)
1340                 {
1341                     /* Replace c endpoint with b */
1342                     swap_em_state(s_b, s_c);
1343                     c   = b;
1344                     gpc = gpb;
1345                 }
1346                 else
1347                 {
1348                     /* Replace a endpoint with b */
1349                     swap_em_state(s_b, s_a);
1350                     a   = b;
1351                     gpa = gpb;
1352                 }
1353
1354                 /*
1355                  * Stop search as soon as we find a value smaller than the endpoints.
1356                  * Never run more than 20 steps, no matter what.
1357                  */
1358                 nminstep++;
1359             }
1360             while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1361                    (nminstep < 20));
1362
1363             if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1364                 nminstep >= 20)
1365             {
1366                 /* OK. We couldn't find a significantly lower energy.
1367                  * If beta==0 this was steepest descent, and then we give up.
1368                  * If not, set beta=0 and restart with steepest descent before quitting.
1369                  */
1370                 if (beta == 0.0)
1371                 {
1372                     /* Converged */
1373                     converged = TRUE;
1374                     break;
1375                 }
1376                 else
1377                 {
1378                     /* Reset memory before giving up */
1379                     beta = 0.0;
1380                     continue;
1381                 }
1382             }
1383
1384             /* Select min energy state of A & C, put the best in B.
1385              */
1386             if (s_c->epot < s_a->epot)
1387             {
1388                 if (debug)
1389                 {
1390                     fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1391                             s_c->epot, s_a->epot);
1392                 }
1393                 swap_em_state(s_b, s_c);
1394                 gpb = gpc;
1395                 b   = c;
1396             }
1397             else
1398             {
1399                 if (debug)
1400                 {
1401                     fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1402                             s_a->epot, s_c->epot);
1403                 }
1404                 swap_em_state(s_b, s_a);
1405                 gpb = gpa;
1406                 b   = a;
1407             }
1408
1409         }
1410         else
1411         {
1412             if (debug)
1413             {
1414                 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1415                         s_c->epot);
1416             }
1417             swap_em_state(s_b, s_c);
1418             gpb = gpc;
1419             b   = c;
1420         }
1421
1422         /* new search direction */
1423         /* beta = 0 means forget all memory and restart with steepest descents. */
1424         if (nstcg && ((step % nstcg) == 0))
1425         {
1426             beta = 0.0;
1427         }
1428         else
1429         {
1430             /* s_min->fnorm cannot be zero, because then we would have converged
1431              * and broken out.
1432              */
1433
1434             /* Polak-Ribiere update.
1435              * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1436              */
1437             beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1438         }
1439         /* Limit beta to prevent oscillations */
1440         if (fabs(beta) > 5.0)
1441         {
1442             beta = 0.0;
1443         }
1444
1445
1446         /* update positions */
1447         swap_em_state(s_min, s_b);
1448         gpa = gpb;
1449
1450         /* Print it if necessary */
1451         if (MASTER(cr))
1452         {
1453             if (bVerbose)
1454             {
1455                 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1456                         step, s_min->epot, s_min->fnorm/sqrt(state_global->natoms),
1457                         s_min->fmax, s_min->a_fmax+1);
1458             }
1459             /* Store the new (lower) energies */
1460             upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1461                        mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1462                        NULL, NULL, vir, pres, NULL, mu_tot, constr);
1463
1464             do_log = do_per_step(step, inputrec->nstlog);
1465             do_ene = do_per_step(step, inputrec->nstenergy);
1466
1467             /* Prepare IMD energy record, if bIMD is TRUE. */
1468             IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1469
1470             if (do_log)
1471             {
1472                 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1473             }
1474             print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1475                        do_log ? fplog : NULL, step, step, eprNORMAL,
1476                        TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1477         }
1478
1479         /* Send energies and positions to the IMD client if bIMD is TRUE. */
1480         if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
1481         {
1482             IMD_send_positions(inputrec->imd);
1483         }
1484
1485         /* Stop when the maximum force lies below tolerance.
1486          * If we have reached machine precision, converged is already set to true.
1487          */
1488         converged = converged || (s_min->fmax < inputrec->em_tol);
1489
1490     } /* End of the loop */
1491
1492     /* IMD cleanup, if bIMD is TRUE. */
1493     IMD_finalize(inputrec->bIMD, inputrec->imd);
1494
1495     if (converged)
1496     {
1497         step--; /* we never took that last step in this case */
1498
1499     }
1500     if (s_min->fmax > inputrec->em_tol)
1501     {
1502         if (MASTER(cr))
1503         {
1504             warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1505             warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1506         }
1507         converged = FALSE;
1508     }
1509
1510     if (MASTER(cr))
1511     {
1512         /* If we printed energy and/or logfile last step (which was the last step)
1513          * we don't have to do it again, but otherwise print the final values.
1514          */
1515         if (!do_log)
1516         {
1517             /* Write final value to log since we didn't do anything the last step */
1518             print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1519         }
1520         if (!do_ene || !do_log)
1521         {
1522             /* Write final energy file entries */
1523             print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1524                        !do_log ? fplog : NULL, step, step, eprNORMAL,
1525                        TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1526         }
1527     }
1528
1529     /* Print some stuff... */
1530     if (MASTER(cr))
1531     {
1532         fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1533     }
1534
1535     /* IMPORTANT!
1536      * For accurate normal mode calculation it is imperative that we
1537      * store the last conformation into the full precision binary trajectory.
1538      *
1539      * However, we should only do it if we did NOT already write this step
1540      * above (which we did if do_x or do_f was true).
1541      */
1542     do_x = !do_per_step(step, inputrec->nstxout);
1543     do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1544
1545     write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1546                   top_global, inputrec, step,
1547                   s_min, state_global, f_global);
1548
1549     fnormn = s_min->fnorm/sqrt(state_global->natoms);
1550
1551     if (MASTER(cr))
1552     {
1553         print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1554                         s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1555         print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1556                         s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1557
1558         fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1559     }
1560
1561     finish_em(cr, outf, walltime_accounting, wcycle);
1562
1563     /* To print the actual number of steps we needed somewhere */
1564     walltime_accounting_set_nsteps_done(walltime_accounting, step);
1565
1566     return 0;
1567 } /* That's all folks */
1568
1569
1570 double do_lbfgs(FILE *fplog, t_commrec *cr,
1571                 int nfile, const t_filenm fnm[],
1572                 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
1573                 int gmx_unused nstglobalcomm,
1574                 gmx_vsite_t *vsite, gmx_constr_t constr,
1575                 int gmx_unused stepout,
1576                 t_inputrec *inputrec,
1577                 gmx_mtop_t *top_global, t_fcdata *fcd,
1578                 t_state *state,
1579                 t_mdatoms *mdatoms,
1580                 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1581                 gmx_edsam_t gmx_unused ed,
1582                 t_forcerec *fr,
1583                 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1584                 gmx_membed_t gmx_unused membed,
1585                 real gmx_unused cpt_period, real gmx_unused max_hours,
1586                 const char gmx_unused *deviceOptions,
1587                 int imdport,
1588                 unsigned long gmx_unused Flags,
1589                 gmx_walltime_accounting_t walltime_accounting)
1590 {
1591     static const char *LBFGS = "Low-Memory BFGS Minimizer";
1592     em_state_t         ems;
1593     gmx_localtop_t    *top;
1594     gmx_enerdata_t    *enerd;
1595     rvec              *f;
1596     gmx_global_stat_t  gstat;
1597     t_graph           *graph;
1598     rvec              *f_global;
1599     int                ncorr, nmaxcorr, point, cp, neval, nminstep;
1600     double             stepsize, gpa, gpb, gpc, tmp, minstep;
1601     real              *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1602     real              *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1603     real               a, b, c, maxdelta, delta;
1604     real               diag, Epot0, Epot, EpotA, EpotB, EpotC;
1605     real               dgdx, dgdg, sq, yr, beta;
1606     t_mdebin          *mdebin;
1607     gmx_bool           converged, first;
1608     rvec               mu_tot;
1609     real               fnorm, fmax;
1610     gmx_bool           do_log, do_ene, do_x, do_f, foundlower, *frozen;
1611     tensor             vir, pres;
1612     int                start, end, number_steps;
1613     gmx_mdoutf_t       outf;
1614     int                i, k, m, n, nfmax, gf, step;
1615     int                mdof_flags;
1616     /* not used */
1617     real               terminate;
1618
1619     if (PAR(cr))
1620     {
1621         gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1622     }
1623
1624     if (NULL != constr)
1625     {
1626         gmx_fatal(FARGS, "The combination of constraints and L-BFGS minimization is not implemented. Either do not use constraints, or use another minimizer (e.g. steepest descent).");
1627     }
1628
1629     n        = 3*state->natoms;
1630     nmaxcorr = inputrec->nbfgscorr;
1631
1632     /* Allocate memory */
1633     /* Use pointers to real so we dont have to loop over both atoms and
1634      * dimensions all the time...
1635      * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1636      * that point to the same memory.
1637      */
1638     snew(xa, n);
1639     snew(xb, n);
1640     snew(xc, n);
1641     snew(fa, n);
1642     snew(fb, n);
1643     snew(fc, n);
1644     snew(frozen, n);
1645
1646     snew(p, n);
1647     snew(lastx, n);
1648     snew(lastf, n);
1649     snew(rho, nmaxcorr);
1650     snew(alpha, nmaxcorr);
1651
1652     snew(dx, nmaxcorr);
1653     for (i = 0; i < nmaxcorr; i++)
1654     {
1655         snew(dx[i], n);
1656     }
1657
1658     snew(dg, nmaxcorr);
1659     for (i = 0; i < nmaxcorr; i++)
1660     {
1661         snew(dg[i], n);
1662     }
1663
1664     step  = 0;
1665     neval = 0;
1666
1667     /* Init em */
1668     init_em(fplog, LBFGS, cr, inputrec,
1669             state, top_global, &ems, &top, &f, &f_global,
1670             nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1671             nfile, fnm, &outf, &mdebin, imdport, Flags);
1672     /* Do_lbfgs is not completely updated like do_steep and do_cg,
1673      * so we free some memory again.
1674      */
1675     sfree(ems.s.x);
1676     sfree(ems.f);
1677
1678     xx = (real *)state->x;
1679     ff = (real *)f;
1680
1681     start = 0;
1682     end   = mdatoms->homenr;
1683
1684     /* Print to log file */
1685     print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1686
1687     do_log = do_ene = do_x = do_f = TRUE;
1688
1689     /* Max number of steps */
1690     number_steps = inputrec->nsteps;
1691
1692     /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1693     gf = 0;
1694     for (i = start; i < end; i++)
1695     {
1696         if (mdatoms->cFREEZE)
1697         {
1698             gf = mdatoms->cFREEZE[i];
1699         }
1700         for (m = 0; m < DIM; m++)
1701         {
1702             frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1703         }
1704     }
1705     if (MASTER(cr))
1706     {
1707         sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1708     }
1709     if (fplog)
1710     {
1711         sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1712     }
1713
1714     if (vsite)
1715     {
1716         construct_vsites(vsite, state->x, 1, NULL,
1717                          top->idef.iparams, top->idef.il,
1718                          fr->ePBC, fr->bMolPBC, cr, state->box);
1719     }
1720
1721     /* Call the force routine and some auxiliary (neighboursearching etc.) */
1722     /* do_force always puts the charge groups in the box and shifts again
1723      * We do not unshift, so molecules are always whole
1724      */
1725     neval++;
1726     ems.s.x = state->x;
1727     ems.f   = f;
1728     evaluate_energy(fplog, cr,
1729                     top_global, &ems, top,
1730                     inputrec, nrnb, wcycle, gstat,
1731                     vsite, constr, fcd, graph, mdatoms, fr,
1732                     mu_tot, enerd, vir, pres, -1, TRUE);
1733     where();
1734
1735     if (MASTER(cr))
1736     {
1737         /* Copy stuff to the energy bin for easy printing etc. */
1738         upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1739                    mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
1740                    NULL, NULL, vir, pres, NULL, mu_tot, constr);
1741
1742         print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
1743         print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1744                    TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1745     }
1746     where();
1747
1748     /* This is the starting energy */
1749     Epot = enerd->term[F_EPOT];
1750
1751     fnorm = ems.fnorm;
1752     fmax  = ems.fmax;
1753     nfmax = ems.a_fmax;
1754
1755     /* Set the initial step.
1756      * since it will be multiplied by the non-normalized search direction
1757      * vector (force vector the first time), we scale it by the
1758      * norm of the force.
1759      */
1760
1761     if (MASTER(cr))
1762     {
1763         fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1764         fprintf(stderr, "   F-max             = %12.5e on atom %d\n", fmax, nfmax+1);
1765         fprintf(stderr, "   F-Norm            = %12.5e\n", fnorm/sqrt(state->natoms));
1766         fprintf(stderr, "\n");
1767         /* and copy to the log file too... */
1768         fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1769         fprintf(fplog, "   F-max             = %12.5e on atom %d\n", fmax, nfmax+1);
1770         fprintf(fplog, "   F-Norm            = %12.5e\n", fnorm/sqrt(state->natoms));
1771         fprintf(fplog, "\n");
1772     }
1773
1774     point = 0;
1775     for (i = 0; i < n; i++)
1776     {
1777         if (!frozen[i])
1778         {
1779             dx[point][i] = ff[i]; /* Initial search direction */
1780         }
1781         else
1782         {
1783             dx[point][i] = 0;
1784         }
1785     }
1786
1787     stepsize  = 1.0/fnorm;
1788
1789     /* Start the loop over BFGS steps.
1790      * Each successful step is counted, and we continue until
1791      * we either converge or reach the max number of steps.
1792      */
1793
1794     ncorr = 0;
1795
1796     /* Set the gradient from the force */
1797     converged = FALSE;
1798     for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1799     {
1800
1801         /* Write coordinates if necessary */
1802         do_x = do_per_step(step, inputrec->nstxout);
1803         do_f = do_per_step(step, inputrec->nstfout);
1804
1805         mdof_flags = 0;
1806         if (do_x)
1807         {
1808             mdof_flags |= MDOF_X;
1809         }
1810
1811         if (do_f)
1812         {
1813             mdof_flags |= MDOF_F;
1814         }
1815
1816         if (inputrec->bIMD)
1817         {
1818             mdof_flags |= MDOF_IMD;
1819         }
1820
1821         mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1822                                          top_global, step, (real)step, state, state, f, f);
1823
1824         /* Do the linesearching in the direction dx[point][0..(n-1)] */
1825
1826         /* pointer to current direction - point=0 first time here */
1827         s = dx[point];
1828
1829         /* calculate line gradient */
1830         for (gpa = 0, i = 0; i < n; i++)
1831         {
1832             gpa -= s[i]*ff[i];
1833         }
1834
1835         /* Calculate minimum allowed stepsize, before the average (norm)
1836          * relative change in coordinate is smaller than precision
1837          */
1838         for (minstep = 0, i = 0; i < n; i++)
1839         {
1840             tmp = fabs(xx[i]);
1841             if (tmp < 1.0)
1842             {
1843                 tmp = 1.0;
1844             }
1845             tmp      = s[i]/tmp;
1846             minstep += tmp*tmp;
1847         }
1848         minstep = GMX_REAL_EPS/sqrt(minstep/n);
1849
1850         if (stepsize < minstep)
1851         {
1852             converged = TRUE;
1853             break;
1854         }
1855
1856         /* Store old forces and coordinates */
1857         for (i = 0; i < n; i++)
1858         {
1859             lastx[i] = xx[i];
1860             lastf[i] = ff[i];
1861         }
1862         Epot0 = Epot;
1863
1864         first = TRUE;
1865
1866         for (i = 0; i < n; i++)
1867         {
1868             xa[i] = xx[i];
1869         }
1870
1871         /* Take a step downhill.
1872          * In theory, we should minimize the function along this direction.
1873          * That is quite possible, but it turns out to take 5-10 function evaluations
1874          * for each line. However, we dont really need to find the exact minimum -
1875          * it is much better to start a new BFGS step in a modified direction as soon
1876          * as we are close to it. This will save a lot of energy evaluations.
1877          *
1878          * In practice, we just try to take a single step.
1879          * If it worked (i.e. lowered the energy), we increase the stepsize but
1880          * the continue straight to the next BFGS step without trying to find any minimum.
1881          * If it didn't work (higher energy), there must be a minimum somewhere between
1882          * the old position and the new one.
1883          *
1884          * Due to the finite numerical accuracy, it turns out that it is a good idea
1885          * to even accept a SMALL increase in energy, if the derivative is still downhill.
1886          * This leads to lower final energies in the tests I've done. / Erik
1887          */
1888         foundlower = FALSE;
1889         EpotA      = Epot0;
1890         a          = 0.0;
1891         c          = a + stepsize; /* reference position along line is zero */
1892
1893         /* Check stepsize first. We do not allow displacements
1894          * larger than emstep.
1895          */
1896         do
1897         {
1898             c        = a + stepsize;
1899             maxdelta = 0;
1900             for (i = 0; i < n; i++)
1901             {
1902                 delta = c*s[i];
1903                 if (delta > maxdelta)
1904                 {
1905                     maxdelta = delta;
1906                 }
1907             }
1908             if (maxdelta > inputrec->em_stepsize)
1909             {
1910                 stepsize *= 0.1;
1911             }
1912         }
1913         while (maxdelta > inputrec->em_stepsize);
1914
1915         /* Take a trial step */
1916         for (i = 0; i < n; i++)
1917         {
1918             xc[i] = lastx[i] + c*s[i];
1919         }
1920
1921         neval++;
1922         /* Calculate energy for the trial step */
1923         ems.s.x = (rvec *)xc;
1924         ems.f   = (rvec *)fc;
1925         evaluate_energy(fplog, cr,
1926                         top_global, &ems, top,
1927                         inputrec, nrnb, wcycle, gstat,
1928                         vsite, constr, fcd, graph, mdatoms, fr,
1929                         mu_tot, enerd, vir, pres, step, FALSE);
1930         EpotC = ems.epot;
1931
1932         /* Calc derivative along line */
1933         for (gpc = 0, i = 0; i < n; i++)
1934         {
1935             gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1936         }
1937         /* Sum the gradient along the line across CPUs */
1938         if (PAR(cr))
1939         {
1940             gmx_sumd(1, &gpc, cr);
1941         }
1942
1943         /* This is the max amount of increase in energy we tolerate */
1944         tmp = sqrt(GMX_REAL_EPS)*fabs(EpotA);
1945
1946         /* Accept the step if the energy is lower, or if it is not significantly higher
1947          * and the line derivative is still negative.
1948          */
1949         if (EpotC < EpotA || (gpc < 0 && EpotC < (EpotA+tmp)))
1950         {
1951             foundlower = TRUE;
1952             /* Great, we found a better energy. Increase step for next iteration
1953              * if we are still going down, decrease it otherwise
1954              */
1955             if (gpc < 0)
1956             {
1957                 stepsize *= 1.618034; /* The golden section */
1958             }
1959             else
1960             {
1961                 stepsize *= 0.618034; /* 1/golden section */
1962             }
1963         }
1964         else
1965         {
1966             /* New energy is the same or higher. We will have to do some work
1967              * to find a smaller value in the interval. Take smaller step next time!
1968              */
1969             foundlower = FALSE;
1970             stepsize  *= 0.618034;
1971         }
1972
1973         /* OK, if we didn't find a lower value we will have to locate one now - there must
1974          * be one in the interval [a=0,c].
1975          * The same thing is valid here, though: Don't spend dozens of iterations to find
1976          * the line minimum. We try to interpolate based on the derivative at the endpoints,
1977          * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1978          *
1979          * I also have a safeguard for potentially really patological functions so we never
1980          * take more than 20 steps before we give up ...
1981          *
1982          * If we already found a lower value we just skip this step and continue to the update.
1983          */
1984
1985         if (!foundlower)
1986         {
1987
1988             nminstep = 0;
1989             do
1990             {
1991                 /* Select a new trial point.
1992                  * If the derivatives at points a & c have different sign we interpolate to zero,
1993                  * otherwise just do a bisection.
1994                  */
1995
1996                 if (gpa < 0 && gpc > 0)
1997                 {
1998                     b = a + gpa*(a-c)/(gpc-gpa);
1999                 }
2000                 else
2001                 {
2002                     b = 0.5*(a+c);
2003                 }
2004
2005                 /* safeguard if interpolation close to machine accuracy causes errors:
2006                  * never go outside the interval
2007                  */
2008                 if (b <= a || b >= c)
2009                 {
2010                     b = 0.5*(a+c);
2011                 }
2012
2013                 /* Take a trial step */
2014                 for (i = 0; i < n; i++)
2015                 {
2016                     xb[i] = lastx[i] + b*s[i];
2017                 }
2018
2019                 neval++;
2020                 /* Calculate energy for the trial step */
2021                 ems.s.x = (rvec *)xb;
2022                 ems.f   = (rvec *)fb;
2023                 evaluate_energy(fplog, cr,
2024                                 top_global, &ems, top,
2025                                 inputrec, nrnb, wcycle, gstat,
2026                                 vsite, constr, fcd, graph, mdatoms, fr,
2027                                 mu_tot, enerd, vir, pres, step, FALSE);
2028                 EpotB = ems.epot;
2029
2030                 fnorm = ems.fnorm;
2031
2032                 for (gpb = 0, i = 0; i < n; i++)
2033                 {
2034                     gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2035
2036                 }
2037                 /* Sum the gradient along the line across CPUs */
2038                 if (PAR(cr))
2039                 {
2040                     gmx_sumd(1, &gpb, cr);
2041                 }
2042
2043                 /* Keep one of the intervals based on the value of the derivative at the new point */
2044                 if (gpb > 0)
2045                 {
2046                     /* Replace c endpoint with b */
2047                     EpotC = EpotB;
2048                     c     = b;
2049                     gpc   = gpb;
2050                     /* swap coord pointers b/c */
2051                     xtmp = xb;
2052                     ftmp = fb;
2053                     xb   = xc;
2054                     fb   = fc;
2055                     xc   = xtmp;
2056                     fc   = ftmp;
2057                 }
2058                 else
2059                 {
2060                     /* Replace a endpoint with b */
2061                     EpotA = EpotB;
2062                     a     = b;
2063                     gpa   = gpb;
2064                     /* swap coord pointers a/b */
2065                     xtmp = xb;
2066                     ftmp = fb;
2067                     xb   = xa;
2068                     fb   = fa;
2069                     xa   = xtmp;
2070                     fa   = ftmp;
2071                 }
2072
2073                 /*
2074                  * Stop search as soon as we find a value smaller than the endpoints,
2075                  * or if the tolerance is below machine precision.
2076                  * Never run more than 20 steps, no matter what.
2077                  */
2078                 nminstep++;
2079             }
2080             while ((EpotB > EpotA || EpotB > EpotC) && (nminstep < 20));
2081
2082             if (fabs(EpotB-Epot0) < GMX_REAL_EPS || nminstep >= 20)
2083             {
2084                 /* OK. We couldn't find a significantly lower energy.
2085                  * If ncorr==0 this was steepest descent, and then we give up.
2086                  * If not, reset memory to restart as steepest descent before quitting.
2087                  */
2088                 if (ncorr == 0)
2089                 {
2090                     /* Converged */
2091                     converged = TRUE;
2092                     break;
2093                 }
2094                 else
2095                 {
2096                     /* Reset memory */
2097                     ncorr = 0;
2098                     /* Search in gradient direction */
2099                     for (i = 0; i < n; i++)
2100                     {
2101                         dx[point][i] = ff[i];
2102                     }
2103                     /* Reset stepsize */
2104                     stepsize = 1.0/fnorm;
2105                     continue;
2106                 }
2107             }
2108
2109             /* Select min energy state of A & C, put the best in xx/ff/Epot
2110              */
2111             if (EpotC < EpotA)
2112             {
2113                 Epot = EpotC;
2114                 /* Use state C */
2115                 for (i = 0; i < n; i++)
2116                 {
2117                     xx[i] = xc[i];
2118                     ff[i] = fc[i];
2119                 }
2120                 stepsize = c;
2121             }
2122             else
2123             {
2124                 Epot = EpotA;
2125                 /* Use state A */
2126                 for (i = 0; i < n; i++)
2127                 {
2128                     xx[i] = xa[i];
2129                     ff[i] = fa[i];
2130                 }
2131                 stepsize = a;
2132             }
2133
2134         }
2135         else
2136         {
2137             /* found lower */
2138             Epot = EpotC;
2139             /* Use state C */
2140             for (i = 0; i < n; i++)
2141             {
2142                 xx[i] = xc[i];
2143                 ff[i] = fc[i];
2144             }
2145             stepsize = c;
2146         }
2147
2148         /* Update the memory information, and calculate a new
2149          * approximation of the inverse hessian
2150          */
2151
2152         /* Have new data in Epot, xx, ff */
2153         if (ncorr < nmaxcorr)
2154         {
2155             ncorr++;
2156         }
2157
2158         for (i = 0; i < n; i++)
2159         {
2160             dg[point][i]  = lastf[i]-ff[i];
2161             dx[point][i] *= stepsize;
2162         }
2163
2164         dgdg = 0;
2165         dgdx = 0;
2166         for (i = 0; i < n; i++)
2167         {
2168             dgdg += dg[point][i]*dg[point][i];
2169             dgdx += dg[point][i]*dx[point][i];
2170         }
2171
2172         diag = dgdx/dgdg;
2173
2174         rho[point] = 1.0/dgdx;
2175         point++;
2176
2177         if (point >= nmaxcorr)
2178         {
2179             point = 0;
2180         }
2181
2182         /* Update */
2183         for (i = 0; i < n; i++)
2184         {
2185             p[i] = ff[i];
2186         }
2187
2188         cp = point;
2189
2190         /* Recursive update. First go back over the memory points */
2191         for (k = 0; k < ncorr; k++)
2192         {
2193             cp--;
2194             if (cp < 0)
2195             {
2196                 cp = ncorr-1;
2197             }
2198
2199             sq = 0;
2200             for (i = 0; i < n; i++)
2201             {
2202                 sq += dx[cp][i]*p[i];
2203             }
2204
2205             alpha[cp] = rho[cp]*sq;
2206
2207             for (i = 0; i < n; i++)
2208             {
2209                 p[i] -= alpha[cp]*dg[cp][i];
2210             }
2211         }
2212
2213         for (i = 0; i < n; i++)
2214         {
2215             p[i] *= diag;
2216         }
2217
2218         /* And then go forward again */
2219         for (k = 0; k < ncorr; k++)
2220         {
2221             yr = 0;
2222             for (i = 0; i < n; i++)
2223             {
2224                 yr += p[i]*dg[cp][i];
2225             }
2226
2227             beta = rho[cp]*yr;
2228             beta = alpha[cp]-beta;
2229
2230             for (i = 0; i < n; i++)
2231             {
2232                 p[i] += beta*dx[cp][i];
2233             }
2234
2235             cp++;
2236             if (cp >= ncorr)
2237             {
2238                 cp = 0;
2239             }
2240         }
2241
2242         for (i = 0; i < n; i++)
2243         {
2244             if (!frozen[i])
2245             {
2246                 dx[point][i] = p[i];
2247             }
2248             else
2249             {
2250                 dx[point][i] = 0;
2251             }
2252         }
2253
2254         stepsize = 1.0;
2255
2256         /* Test whether the convergence criterion is met */
2257         get_f_norm_max(cr, &(inputrec->opts), mdatoms, f, &fnorm, &fmax, &nfmax);
2258
2259         /* Print it if necessary */
2260         if (MASTER(cr))
2261         {
2262             if (bVerbose)
2263             {
2264                 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2265                         step, Epot, fnorm/sqrt(state->natoms), fmax, nfmax+1);
2266             }
2267             /* Store the new (lower) energies */
2268             upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2269                        mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
2270                        NULL, NULL, vir, pres, NULL, mu_tot, constr);
2271             do_log = do_per_step(step, inputrec->nstlog);
2272             do_ene = do_per_step(step, inputrec->nstenergy);
2273             if (do_log)
2274             {
2275                 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2276             }
2277             print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2278                        do_log ? fplog : NULL, step, step, eprNORMAL,
2279                        TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2280         }
2281
2282         /* Send x and E to IMD client, if bIMD is TRUE. */
2283         if (do_IMD(inputrec->bIMD, step, cr, TRUE, state->box, state->x, inputrec, 0, wcycle) && MASTER(cr))
2284         {
2285             IMD_send_positions(inputrec->imd);
2286         }
2287
2288         /* Stop when the maximum force lies below tolerance.
2289          * If we have reached machine precision, converged is already set to true.
2290          */
2291
2292         converged = converged || (fmax < inputrec->em_tol);
2293
2294     } /* End of the loop */
2295
2296     /* IMD cleanup, if bIMD is TRUE. */
2297     IMD_finalize(inputrec->bIMD, inputrec->imd);
2298
2299     if (converged)
2300     {
2301         step--; /* we never took that last step in this case */
2302
2303     }
2304     if (fmax > inputrec->em_tol)
2305     {
2306         if (MASTER(cr))
2307         {
2308             warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2309             warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2310         }
2311         converged = FALSE;
2312     }
2313
2314     /* If we printed energy and/or logfile last step (which was the last step)
2315      * we don't have to do it again, but otherwise print the final values.
2316      */
2317     if (!do_log) /* Write final value to log since we didn't do anythin last step */
2318     {
2319         print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2320     }
2321     if (!do_ene || !do_log) /* Write final energy file entries */
2322     {
2323         print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2324                    !do_log ? fplog : NULL, step, step, eprNORMAL,
2325                    TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2326     }
2327
2328     /* Print some stuff... */
2329     if (MASTER(cr))
2330     {
2331         fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2332     }
2333
2334     /* IMPORTANT!
2335      * For accurate normal mode calculation it is imperative that we
2336      * store the last conformation into the full precision binary trajectory.
2337      *
2338      * However, we should only do it if we did NOT already write this step
2339      * above (which we did if do_x or do_f was true).
2340      */
2341     do_x = !do_per_step(step, inputrec->nstxout);
2342     do_f = !do_per_step(step, inputrec->nstfout);
2343     write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2344                   top_global, inputrec, step,
2345                   &ems, state, f);
2346
2347     if (MASTER(cr))
2348     {
2349         print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2350                         number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2351         print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2352                         number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2353
2354         fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2355     }
2356
2357     finish_em(cr, outf, walltime_accounting, wcycle);
2358
2359     /* To print the actual number of steps we needed somewhere */
2360     walltime_accounting_set_nsteps_done(walltime_accounting, step);
2361
2362     return 0;
2363 } /* That's all folks */
2364
2365
2366 double do_steep(FILE *fplog, t_commrec *cr,
2367                 int nfile, const t_filenm fnm[],
2368                 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2369                 int gmx_unused nstglobalcomm,
2370                 gmx_vsite_t *vsite, gmx_constr_t constr,
2371                 int gmx_unused stepout,
2372                 t_inputrec *inputrec,
2373                 gmx_mtop_t *top_global, t_fcdata *fcd,
2374                 t_state *state_global,
2375                 t_mdatoms *mdatoms,
2376                 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2377                 gmx_edsam_t gmx_unused  ed,
2378                 t_forcerec *fr,
2379                 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2380                 gmx_membed_t gmx_unused membed,
2381                 real gmx_unused cpt_period, real gmx_unused max_hours,
2382                 const char  gmx_unused *deviceOptions,
2383                 int imdport,
2384                 unsigned long gmx_unused Flags,
2385                 gmx_walltime_accounting_t walltime_accounting)
2386 {
2387     const char       *SD = "Steepest Descents";
2388     em_state_t       *s_min, *s_try;
2389     rvec             *f_global;
2390     gmx_localtop_t   *top;
2391     gmx_enerdata_t   *enerd;
2392     rvec             *f;
2393     gmx_global_stat_t gstat;
2394     t_graph          *graph;
2395     real              stepsize, constepsize;
2396     real              ustep, fnormn;
2397     gmx_mdoutf_t      outf;
2398     t_mdebin         *mdebin;
2399     gmx_bool          bDone, bAbort, do_x, do_f;
2400     tensor            vir, pres;
2401     rvec              mu_tot;
2402     int               nsteps;
2403     int               count          = 0;
2404     int               steps_accepted = 0;
2405     /* not used */
2406     real              terminate = 0;
2407
2408     s_min = init_em_state();
2409     s_try = init_em_state();
2410
2411     /* Init em and store the local state in s_try */
2412     init_em(fplog, SD, cr, inputrec,
2413             state_global, top_global, s_try, &top, &f, &f_global,
2414             nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2415             nfile, fnm, &outf, &mdebin, imdport, Flags);
2416
2417     /* Print to log file  */
2418     print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2419
2420     /* Set variables for stepsize (in nm). This is the largest
2421      * step that we are going to make in any direction.
2422      */
2423     ustep    = inputrec->em_stepsize;
2424     stepsize = 0;
2425
2426     /* Max number of steps  */
2427     nsteps = inputrec->nsteps;
2428
2429     if (MASTER(cr))
2430     {
2431         /* Print to the screen  */
2432         sp_header(stderr, SD, inputrec->em_tol, nsteps);
2433     }
2434     if (fplog)
2435     {
2436         sp_header(fplog, SD, inputrec->em_tol, nsteps);
2437     }
2438
2439     /**** HERE STARTS THE LOOP ****
2440      * count is the counter for the number of steps
2441      * bDone will be TRUE when the minimization has converged
2442      * bAbort will be TRUE when nsteps steps have been performed or when
2443      * the stepsize becomes smaller than is reasonable for machine precision
2444      */
2445     count  = 0;
2446     bDone  = FALSE;
2447     bAbort = FALSE;
2448     while (!bDone && !bAbort)
2449     {
2450         bAbort = (nsteps >= 0) && (count == nsteps);
2451
2452         /* set new coordinates, except for first step */
2453         if (count > 0)
2454         {
2455             do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2456                        s_min, stepsize, s_min->f, s_try,
2457                        constr, top, nrnb, wcycle, count);
2458         }
2459
2460         evaluate_energy(fplog, cr,
2461                         top_global, s_try, top,
2462                         inputrec, nrnb, wcycle, gstat,
2463                         vsite, constr, fcd, graph, mdatoms, fr,
2464                         mu_tot, enerd, vir, pres, count, count == 0);
2465
2466         if (MASTER(cr))
2467         {
2468             print_ebin_header(fplog, count, count, s_try->s.lambda[efptFEP]);
2469         }
2470
2471         if (count == 0)
2472         {
2473             s_min->epot = s_try->epot + 1;
2474         }
2475
2476         /* Print it if necessary  */
2477         if (MASTER(cr))
2478         {
2479             if (bVerbose)
2480             {
2481                 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2482                         count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2483                         (s_try->epot < s_min->epot) ? '\n' : '\r');
2484             }
2485
2486             if (s_try->epot < s_min->epot)
2487             {
2488                 /* Store the new (lower) energies  */
2489                 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2490                            mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2491                            s_try->s.box, NULL, NULL, vir, pres, NULL, mu_tot, constr);
2492
2493                 /* Prepare IMD energy record, if bIMD is TRUE. */
2494                 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2495
2496                 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2497                            do_per_step(steps_accepted, inputrec->nstdisreout),
2498                            do_per_step(steps_accepted, inputrec->nstorireout),
2499                            fplog, count, count, eprNORMAL, TRUE,
2500                            mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2501                 fflush(fplog);
2502             }
2503         }
2504
2505         /* Now if the new energy is smaller than the previous...
2506          * or if this is the first step!
2507          * or if we did random steps!
2508          */
2509
2510         if ( (count == 0) || (s_try->epot < s_min->epot) )
2511         {
2512             steps_accepted++;
2513
2514             /* Test whether the convergence criterion is met...  */
2515             bDone = (s_try->fmax < inputrec->em_tol);
2516
2517             /* Copy the arrays for force, positions and energy  */
2518             /* The 'Min' array always holds the coords and forces of the minimal
2519                sampled energy  */
2520             swap_em_state(s_min, s_try);
2521             if (count > 0)
2522             {
2523                 ustep *= 1.2;
2524             }
2525
2526             /* Write to trn, if necessary */
2527             do_x = do_per_step(steps_accepted, inputrec->nstxout);
2528             do_f = do_per_step(steps_accepted, inputrec->nstfout);
2529             write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
2530                           top_global, inputrec, count,
2531                           s_min, state_global, f_global);
2532         }
2533         else
2534         {
2535             /* If energy is not smaller make the step smaller...  */
2536             ustep *= 0.5;
2537
2538             if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2539             {
2540                 /* Reload the old state */
2541                 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2542                                        s_min, top, mdatoms, fr, vsite, constr,
2543                                        nrnb, wcycle);
2544             }
2545         }
2546
2547         /* Determine new step  */
2548         stepsize = ustep/s_min->fmax;
2549
2550         /* Check if stepsize is too small, with 1 nm as a characteristic length */
2551 #ifdef GMX_DOUBLE
2552         if (count == nsteps || ustep < 1e-12)
2553 #else
2554         if (count == nsteps || ustep < 1e-6)
2555 #endif
2556         {
2557             if (MASTER(cr))
2558             {
2559                 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2560                 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2561             }
2562             bAbort = TRUE;
2563         }
2564
2565         /* Send IMD energies and positions, if bIMD is TRUE. */
2566         if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
2567         {
2568             IMD_send_positions(inputrec->imd);
2569         }
2570
2571         count++;
2572     } /* End of the loop  */
2573
2574     /* IMD cleanup, if bIMD is TRUE. */
2575     IMD_finalize(inputrec->bIMD, inputrec->imd);
2576
2577     /* Print some data...  */
2578     if (MASTER(cr))
2579     {
2580         fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2581     }
2582     write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2583                   top_global, inputrec, count,
2584                   s_min, state_global, f_global);
2585
2586     fnormn = s_min->fnorm/sqrt(state_global->natoms);
2587
2588     if (MASTER(cr))
2589     {
2590         print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2591                         s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2592         print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2593                         s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2594     }
2595
2596     finish_em(cr, outf, walltime_accounting, wcycle);
2597
2598     /* To print the actual number of steps we needed somewhere */
2599     inputrec->nsteps = count;
2600
2601     walltime_accounting_set_nsteps_done(walltime_accounting, count);
2602
2603     return 0;
2604 } /* That's all folks */
2605
2606
2607 double do_nm(FILE *fplog, t_commrec *cr,
2608              int nfile, const t_filenm fnm[],
2609              const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused  bCompact,
2610              int gmx_unused nstglobalcomm,
2611              gmx_vsite_t *vsite, gmx_constr_t constr,
2612              int gmx_unused stepout,
2613              t_inputrec *inputrec,
2614              gmx_mtop_t *top_global, t_fcdata *fcd,
2615              t_state *state_global,
2616              t_mdatoms *mdatoms,
2617              t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2618              gmx_edsam_t  gmx_unused ed,
2619              t_forcerec *fr,
2620              int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2621              gmx_membed_t gmx_unused membed,
2622              real gmx_unused cpt_period, real gmx_unused max_hours,
2623              const char gmx_unused *deviceOptions,
2624              int imdport,
2625              unsigned long gmx_unused Flags,
2626              gmx_walltime_accounting_t walltime_accounting)
2627 {
2628     const char          *NM = "Normal Mode Analysis";
2629     gmx_mdoutf_t         outf;
2630     int                  natoms, atom, d;
2631     int                  nnodes, node;
2632     rvec                *f_global;
2633     gmx_localtop_t      *top;
2634     gmx_enerdata_t      *enerd;
2635     rvec                *f;
2636     gmx_global_stat_t    gstat;
2637     t_graph             *graph;
2638     real                 t, t0, lambda, lam0;
2639     gmx_bool             bNS;
2640     tensor               vir, pres;
2641     rvec                 mu_tot;
2642     rvec                *fneg, *dfdx;
2643     gmx_bool             bSparse; /* use sparse matrix storage format */
2644     size_t               sz = 0;
2645     gmx_sparsematrix_t * sparse_matrix           = NULL;
2646     real           *     full_matrix             = NULL;
2647     em_state_t       *   state_work;
2648
2649     /* added with respect to mdrun */
2650     int        i, j, k, row, col;
2651     real       der_range = 10.0*sqrt(GMX_REAL_EPS);
2652     real       x_min;
2653     real       fnorm, fmax;
2654
2655     if (constr != NULL)
2656     {
2657         gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2658     }
2659
2660     state_work = init_em_state();
2661
2662     /* Init em and store the local state in state_minimum */
2663     init_em(fplog, NM, cr, inputrec,
2664             state_global, top_global, state_work, &top,
2665             &f, &f_global,
2666             nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2667             nfile, fnm, &outf, NULL, imdport, Flags);
2668
2669     natoms = top_global->natoms;
2670     snew(fneg, natoms);
2671     snew(dfdx, natoms);
2672
2673 #ifndef GMX_DOUBLE
2674     if (MASTER(cr))
2675     {
2676         fprintf(stderr,
2677                 "NOTE: This version of Gromacs has been compiled in single precision,\n"
2678                 "      which MIGHT not be accurate enough for normal mode analysis.\n"
2679                 "      Gromacs now uses sparse matrix storage, so the memory requirements\n"
2680                 "      are fairly modest even if you recompile in double precision.\n\n");
2681     }
2682 #endif
2683
2684     /* Check if we can/should use sparse storage format.
2685      *
2686      * Sparse format is only useful when the Hessian itself is sparse, which it
2687      * will be when we use a cutoff.
2688      * For small systems (n<1000) it is easier to always use full matrix format, though.
2689      */
2690     if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2691     {
2692         md_print_info(cr, fplog, "Non-cutoff electrostatics used, forcing full Hessian format.\n");
2693         bSparse = FALSE;
2694     }
2695     else if (top_global->natoms < 1000)
2696     {
2697         md_print_info(cr, fplog, "Small system size (N=%d), using full Hessian format.\n", top_global->natoms);
2698         bSparse = FALSE;
2699     }
2700     else
2701     {
2702         md_print_info(cr, fplog, "Using compressed symmetric sparse Hessian format.\n");
2703         bSparse = TRUE;
2704     }
2705
2706     if (MASTER(cr))
2707     {
2708         sz = DIM*top_global->natoms;
2709
2710         fprintf(stderr, "Allocating Hessian memory...\n\n");
2711
2712         if (bSparse)
2713         {
2714             sparse_matrix = gmx_sparsematrix_init(sz);
2715             sparse_matrix->compressed_symmetric = TRUE;
2716         }
2717         else
2718         {
2719             snew(full_matrix, sz*sz);
2720         }
2721     }
2722
2723     /* Initial values */
2724     t0           = inputrec->init_t;
2725     lam0         = inputrec->fepvals->init_lambda;
2726     t            = t0;
2727     lambda       = lam0;
2728
2729     init_nrnb(nrnb);
2730
2731     where();
2732
2733     /* Write start time and temperature */
2734     print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2735
2736     /* fudge nr of steps to nr of atoms */
2737     inputrec->nsteps = natoms*2;
2738
2739     if (MASTER(cr))
2740     {
2741         fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2742                 *(top_global->name), (int)inputrec->nsteps);
2743     }
2744
2745     nnodes = cr->nnodes;
2746
2747     /* Make evaluate_energy do a single node force calculation */
2748     cr->nnodes = 1;
2749     evaluate_energy(fplog, cr,
2750                     top_global, state_work, top,
2751                     inputrec, nrnb, wcycle, gstat,
2752                     vsite, constr, fcd, graph, mdatoms, fr,
2753                     mu_tot, enerd, vir, pres, -1, TRUE);
2754     cr->nnodes = nnodes;
2755
2756     /* if forces are not small, warn user */
2757     get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, state_work);
2758
2759     md_print_info(cr, fplog, "Maximum force:%12.5e\n", state_work->fmax);
2760     if (state_work->fmax > 1.0e-3)
2761     {
2762         md_print_info(cr, fplog,
2763                       "The force is probably not small enough to "
2764                       "ensure that you are at a minimum.\n"
2765                       "Be aware that negative eigenvalues may occur\n"
2766                       "when the resulting matrix is diagonalized.\n\n");
2767     }
2768
2769     /***********************************************************
2770      *
2771      *      Loop over all pairs in matrix
2772      *
2773      *      do_force called twice. Once with positive and
2774      *      once with negative displacement
2775      *
2776      ************************************************************/
2777
2778     /* Steps are divided one by one over the nodes */
2779     for (atom = cr->nodeid; atom < natoms; atom += nnodes)
2780     {
2781
2782         for (d = 0; d < DIM; d++)
2783         {
2784             x_min = state_work->s.x[atom][d];
2785
2786             state_work->s.x[atom][d] = x_min - der_range;
2787
2788             /* Make evaluate_energy do a single node force calculation */
2789             cr->nnodes = 1;
2790             evaluate_energy(fplog, cr,
2791                             top_global, state_work, top,
2792                             inputrec, nrnb, wcycle, gstat,
2793                             vsite, constr, fcd, graph, mdatoms, fr,
2794                             mu_tot, enerd, vir, pres, atom*2, FALSE);
2795
2796             for (i = 0; i < natoms; i++)
2797             {
2798                 copy_rvec(state_work->f[i], fneg[i]);
2799             }
2800
2801             state_work->s.x[atom][d] = x_min + der_range;
2802
2803             evaluate_energy(fplog, cr,
2804                             top_global, state_work, top,
2805                             inputrec, nrnb, wcycle, gstat,
2806                             vsite, constr, fcd, graph, mdatoms, fr,
2807                             mu_tot, enerd, vir, pres, atom*2+1, FALSE);
2808             cr->nnodes = nnodes;
2809
2810             /* x is restored to original */
2811             state_work->s.x[atom][d] = x_min;
2812
2813             for (j = 0; j < natoms; j++)
2814             {
2815                 for (k = 0; (k < DIM); k++)
2816                 {
2817                     dfdx[j][k] =
2818                         -(state_work->f[j][k] - fneg[j][k])/(2*der_range);
2819                 }
2820             }
2821
2822             if (!MASTER(cr))
2823             {
2824 #ifdef GMX_MPI
2825 #ifdef GMX_DOUBLE
2826 #define mpi_type MPI_DOUBLE
2827 #else
2828 #define mpi_type MPI_FLOAT
2829 #endif
2830                 MPI_Send(dfdx[0], natoms*DIM, mpi_type, MASTERNODE(cr), cr->nodeid,
2831                          cr->mpi_comm_mygroup);
2832 #endif
2833             }
2834             else
2835             {
2836                 for (node = 0; (node < nnodes && atom+node < natoms); node++)
2837                 {
2838                     if (node > 0)
2839                     {
2840 #ifdef GMX_MPI
2841                         MPI_Status stat;
2842                         MPI_Recv(dfdx[0], natoms*DIM, mpi_type, node, node,
2843                                  cr->mpi_comm_mygroup, &stat);
2844 #undef mpi_type
2845 #endif
2846                     }
2847
2848                     row = (atom + node)*DIM + d;
2849
2850                     for (j = 0; j < natoms; j++)
2851                     {
2852                         for (k = 0; k < DIM; k++)
2853                         {
2854                             col = j*DIM + k;
2855
2856                             if (bSparse)
2857                             {
2858                                 if (col >= row && dfdx[j][k] != 0.0)
2859                                 {
2860                                     gmx_sparsematrix_increment_value(sparse_matrix,
2861                                                                      row, col, dfdx[j][k]);
2862                                 }
2863                             }
2864                             else
2865                             {
2866                                 full_matrix[row*sz+col] = dfdx[j][k];
2867                             }
2868                         }
2869                     }
2870                 }
2871             }
2872
2873             if (bVerbose && fplog)
2874             {
2875                 fflush(fplog);
2876             }
2877         }
2878         /* write progress */
2879         if (MASTER(cr) && bVerbose)
2880         {
2881             fprintf(stderr, "\rFinished step %d out of %d",
2882                     min(atom+nnodes, natoms), natoms);
2883             fflush(stderr);
2884         }
2885     }
2886
2887     if (MASTER(cr))
2888     {
2889         fprintf(stderr, "\n\nWriting Hessian...\n");
2890         gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2891     }
2892
2893     finish_em(cr, outf, walltime_accounting, wcycle);
2894
2895     walltime_accounting_set_nsteps_done(walltime_accounting, natoms*2);
2896
2897     return 0;
2898 }