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