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