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