2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017, 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.
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.
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.
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.
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.
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.
39 * \brief This file defines integrators for energy minimization
41 * \author Berk Hess <hess@kth.se>
42 * \author Erik Lindahl <erik@kth.se>
43 * \ingroup module_mdlib
58 #include "gromacs/commandline/filenm.h"
59 #include "gromacs/domdec/domdec.h"
60 #include "gromacs/domdec/domdec_struct.h"
61 #include "gromacs/ewald/pme.h"
62 #include "gromacs/fileio/confio.h"
63 #include "gromacs/fileio/mtxio.h"
64 #include "gromacs/gmxlib/network.h"
65 #include "gromacs/gmxlib/nrnb.h"
66 #include "gromacs/imd/imd.h"
67 #include "gromacs/linearalgebra/sparsematrix.h"
68 #include "gromacs/listed-forces/manage-threading.h"
69 #include "gromacs/math/functions.h"
70 #include "gromacs/math/vec.h"
71 #include "gromacs/mdlib/constr.h"
72 #include "gromacs/mdlib/force.h"
73 #include "gromacs/mdlib/forcerec.h"
74 #include "gromacs/mdlib/gmx_omp_nthreads.h"
75 #include "gromacs/mdlib/md_support.h"
76 #include "gromacs/mdlib/mdatoms.h"
77 #include "gromacs/mdlib/mdebin.h"
78 #include "gromacs/mdlib/mdrun.h"
79 #include "gromacs/mdlib/mdsetup.h"
80 #include "gromacs/mdlib/ns.h"
81 #include "gromacs/mdlib/shellfc.h"
82 #include "gromacs/mdlib/sim_util.h"
83 #include "gromacs/mdlib/tgroup.h"
84 #include "gromacs/mdlib/trajectory_writing.h"
85 #include "gromacs/mdlib/update.h"
86 #include "gromacs/mdlib/vsite.h"
87 #include "gromacs/mdtypes/commrec.h"
88 #include "gromacs/mdtypes/inputrec.h"
89 #include "gromacs/mdtypes/md_enums.h"
90 #include "gromacs/mdtypes/state.h"
91 #include "gromacs/pbcutil/mshift.h"
92 #include "gromacs/pbcutil/pbc.h"
93 #include "gromacs/timing/wallcycle.h"
94 #include "gromacs/timing/walltime_accounting.h"
95 #include "gromacs/topology/mtop_util.h"
96 #include "gromacs/topology/topology.h"
97 #include "gromacs/utility/cstringutil.h"
98 #include "gromacs/utility/exceptions.h"
99 #include "gromacs/utility/fatalerror.h"
100 #include "gromacs/utility/logger.h"
101 #include "gromacs/utility/smalloc.h"
103 //! Utility structure for manipulating states during EM
105 //! Copy of the global state
111 //! Norm of the force
119 //! Print the EM starting conditions
120 static void print_em_start(FILE *fplog
,
122 gmx_walltime_accounting_t walltime_accounting
,
123 gmx_wallcycle_t wcycle
,
126 walltime_accounting_start(walltime_accounting
);
127 wallcycle_start(wcycle
, ewcRUN
);
128 print_start(fplog
, cr
, walltime_accounting
, name
);
131 //! Stop counting time for EM
132 static void em_time_end(gmx_walltime_accounting_t walltime_accounting
,
133 gmx_wallcycle_t wcycle
)
135 wallcycle_stop(wcycle
, ewcRUN
);
137 walltime_accounting_end(walltime_accounting
);
140 //! Printing a log file and console header
141 static void sp_header(FILE *out
, const char *minimizer
, real ftol
, int nsteps
)
144 fprintf(out
, "%s:\n", minimizer
);
145 fprintf(out
, " Tolerance (Fmax) = %12.5e\n", ftol
);
146 fprintf(out
, " Number of steps = %12d\n", nsteps
);
149 //! Print warning message
150 static void warn_step(FILE *fp
, real ftol
, gmx_bool bLastStep
, gmx_bool bConstrain
)
156 "\nEnergy minimization reached the maximum number "
157 "of steps before the forces reached the requested "
158 "precision Fmax < %g.\n", ftol
);
163 "\nEnergy minimization has stopped, but the forces have "
164 "not converged to the requested precision Fmax < %g (which "
165 "may not be possible for your system). It stopped "
166 "because the algorithm tried to make a new step whose size "
167 "was too small, or there was no change in the energy since "
168 "last step. Either way, we regard the minimization as "
169 "converged to within the available machine precision, "
170 "given your starting configuration and EM parameters.\n%s%s",
172 sizeof(real
) < sizeof(double) ?
173 "\nDouble precision normally gives you higher accuracy, but "
174 "this is often not needed for preparing to run molecular "
178 "You might need to increase your constraint accuracy, or turn\n"
179 "off constraints altogether (set constraints = none in mdp file)\n" :
182 fputs(wrap_lines(buffer
, 78, 0, FALSE
), fp
);
185 //! Print message about convergence of the EM
186 static void print_converged(FILE *fp
, const char *alg
, real ftol
,
187 gmx_int64_t count
, gmx_bool bDone
, gmx_int64_t nsteps
,
188 const em_state_t
*ems
, double sqrtNumAtoms
)
190 char buf
[STEPSTRSIZE
];
194 fprintf(fp
, "\n%s converged to Fmax < %g in %s steps\n",
195 alg
, ftol
, gmx_step_str(count
, buf
));
197 else if (count
< nsteps
)
199 fprintf(fp
, "\n%s converged to machine precision in %s steps,\n"
200 "but did not reach the requested Fmax < %g.\n",
201 alg
, gmx_step_str(count
, buf
), ftol
);
205 fprintf(fp
, "\n%s did not converge to Fmax < %g in %s steps.\n",
206 alg
, ftol
, gmx_step_str(count
, buf
));
210 fprintf(fp
, "Potential Energy = %21.14e\n", ems
->epot
);
211 fprintf(fp
, "Maximum force = %21.14e on atom %d\n", ems
->fmax
, ems
->a_fmax
+ 1);
212 fprintf(fp
, "Norm of force = %21.14e\n", ems
->fnorm
/sqrtNumAtoms
);
214 fprintf(fp
, "Potential Energy = %14.7e\n", ems
->epot
);
215 fprintf(fp
, "Maximum force = %14.7e on atom %d\n", ems
->fmax
, ems
->a_fmax
+ 1);
216 fprintf(fp
, "Norm of force = %14.7e\n", ems
->fnorm
/sqrtNumAtoms
);
220 //! Compute the norm and max of the force array in parallel
221 static void get_f_norm_max(t_commrec
*cr
,
222 t_grpopts
*opts
, t_mdatoms
*mdatoms
, const rvec
*f
,
223 real
*fnorm
, real
*fmax
, int *a_fmax
)
227 int la_max
, a_max
, start
, end
, i
, m
, gf
;
229 /* This routine finds the largest force and returns it.
230 * On parallel machines the global max is taken.
236 end
= mdatoms
->homenr
;
237 if (mdatoms
->cFREEZE
)
239 for (i
= start
; i
< end
; i
++)
241 gf
= mdatoms
->cFREEZE
[i
];
243 for (m
= 0; m
< DIM
; m
++)
245 if (!opts
->nFreeze
[gf
][m
])
247 fam
+= gmx::square(f
[i
][m
]);
260 for (i
= start
; i
< end
; i
++)
272 if (la_max
>= 0 && DOMAINDECOMP(cr
))
274 a_max
= cr
->dd
->gatindex
[la_max
];
282 snew(sum
, 2*cr
->nnodes
+1);
283 sum
[2*cr
->nodeid
] = fmax2
;
284 sum
[2*cr
->nodeid
+1] = a_max
;
285 sum
[2*cr
->nnodes
] = fnorm2
;
286 gmx_sumd(2*cr
->nnodes
+1, sum
, cr
);
287 fnorm2
= sum
[2*cr
->nnodes
];
288 /* Determine the global maximum */
289 for (i
= 0; i
< cr
->nnodes
; i
++)
291 if (sum
[2*i
] > fmax2
)
294 a_max
= (int)(sum
[2*i
+1] + 0.5);
302 *fnorm
= sqrt(fnorm2
);
314 //! Compute the norm of the force
315 static void get_state_f_norm_max(t_commrec
*cr
,
316 t_grpopts
*opts
, t_mdatoms
*mdatoms
,
319 get_f_norm_max(cr
, opts
, mdatoms
, as_rvec_array(ems
->f
.data()),
320 &ems
->fnorm
, &ems
->fmax
, &ems
->a_fmax
);
323 //! Initialize the energy minimization
324 static void init_em(FILE *fplog
, const char *title
,
325 t_commrec
*cr
, gmx::IMDOutputProvider
*outputProvider
,
327 t_state
*state_global
, gmx_mtop_t
*top_global
,
328 em_state_t
*ems
, gmx_localtop_t
**top
,
329 t_nrnb
*nrnb
, rvec mu_tot
,
330 t_forcerec
*fr
, gmx_enerdata_t
**enerd
,
331 t_graph
**graph
, t_mdatoms
*mdatoms
, gmx_global_stat_t
*gstat
,
332 gmx_vsite_t
*vsite
, gmx_constr_t constr
, gmx_shellfc_t
**shellfc
,
333 int nfile
, const t_filenm fnm
[],
334 gmx_mdoutf_t
*outf
, t_mdebin
**mdebin
,
335 int imdport
, unsigned long gmx_unused Flags
,
336 gmx_wallcycle_t wcycle
)
342 fprintf(fplog
, "Initiating %s\n", title
);
345 state_global
->ngtc
= 0;
347 /* Initialize lambda variables */
348 initialize_lambdas(fplog
, ir
, &(state_global
->fep_state
), state_global
->lambda
, nullptr);
352 /* Interactive molecular dynamics */
353 init_IMD(ir
, cr
, top_global
, fplog
, 1, as_rvec_array(state_global
->x
.data()),
354 nfile
, fnm
, nullptr, imdport
, Flags
);
358 GMX_ASSERT(shellfc
!= NULL
, "With NM we always support shells");
360 *shellfc
= init_shell_flexcon(stdout
,
362 n_flexible_constraints(constr
),
368 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir
->eI
), "This else currently only handles energy minimizers, consider if your algorithm needs shell/flexible-constraint support");
370 /* With energy minimization, shells and flexible constraints are
371 * automatically minimized when treated like normal DOFS.
373 if (shellfc
!= nullptr)
379 if (DOMAINDECOMP(cr
))
381 *top
= dd_init_local_top(top_global
);
383 dd_init_local_state(cr
->dd
, state_global
, &ems
->s
);
385 /* Distribute the charge groups over the nodes from the master node */
386 dd_partition_system(fplog
, ir
->init_step
, cr
, TRUE
, 1,
387 state_global
, top_global
, ir
,
388 &ems
->s
, &ems
->f
, mdatoms
, *top
,
390 nrnb
, nullptr, FALSE
);
391 dd_store_state(cr
->dd
, &ems
->s
);
397 state_change_natoms(state_global
, state_global
->natoms
);
398 /* Just copy the state */
399 ems
->s
= *state_global
;
400 state_change_natoms(&ems
->s
, ems
->s
.natoms
);
401 /* We need to allocate one element extra, since we might use
402 * (unaligned) 4-wide SIMD loads to access rvec entries.
404 ems
->f
.resize(ems
->s
.natoms
+ 1);
407 mdAlgorithmsSetupAtomData(cr
, ir
, top_global
, *top
, fr
,
409 vsite
, shellfc
? *shellfc
: nullptr);
413 set_vsite_top(vsite
, *top
, mdatoms
, cr
);
417 update_mdatoms(mdatoms
, state_global
->lambda
[efptMASS
]);
421 if (ir
->eConstrAlg
== econtSHAKE
&&
422 gmx_mtop_ftype_count(top_global
, F_CONSTR
) > 0)
424 gmx_fatal(FARGS
, "Can not do energy minimization with %s, use %s\n",
425 econstr_names
[econtSHAKE
], econstr_names
[econtLINCS
]);
428 if (!DOMAINDECOMP(cr
))
430 set_constraints(constr
, *top
, ir
, mdatoms
, cr
);
433 if (!ir
->bContinuation
)
435 /* Constrain the starting coordinates */
437 constrain(PAR(cr
) ? nullptr : fplog
, TRUE
, TRUE
, constr
, &(*top
)->idef
,
438 ir
, cr
, -1, 0, 1.0, mdatoms
,
439 as_rvec_array(ems
->s
.x
.data()),
440 as_rvec_array(ems
->s
.x
.data()),
442 fr
->bMolPBC
, ems
->s
.box
,
443 ems
->s
.lambda
[efptFEP
], &dvdl_constr
,
444 nullptr, nullptr, nrnb
, econqCoord
);
450 *gstat
= global_stat_init(ir
);
457 *outf
= init_mdoutf(fplog
, nfile
, fnm
, 0, cr
, outputProvider
, ir
, top_global
, nullptr, wcycle
);
460 init_enerdata(top_global
->groups
.grps
[egcENER
].nr
, ir
->fepvals
->n_lambda
,
463 if (mdebin
!= nullptr)
465 /* Init bin for energy stuff */
466 *mdebin
= init_mdebin(mdoutf_get_fp_ene(*outf
), top_global
, ir
, nullptr);
470 calc_shifts(ems
->s
.box
, fr
->shift_vec
);
473 //! Finalize the minimization
474 static void finish_em(t_commrec
*cr
, gmx_mdoutf_t outf
,
475 gmx_walltime_accounting_t walltime_accounting
,
476 gmx_wallcycle_t wcycle
)
478 if (!(cr
->duty
& DUTY_PME
))
480 /* Tell the PME only node to finish */
481 gmx_pme_send_finish(cr
);
486 em_time_end(walltime_accounting
, wcycle
);
489 //! Swap two different EM states during minimization
490 static void swap_em_state(em_state_t
**ems1
, em_state_t
**ems2
)
499 //! Save the EM trajectory
500 static void write_em_traj(FILE *fplog
, t_commrec
*cr
,
502 gmx_bool bX
, gmx_bool bF
, const char *confout
,
503 gmx_mtop_t
*top_global
,
504 t_inputrec
*ir
, gmx_int64_t step
,
506 t_state
*state_global
,
507 ObservablesHistory
*observablesHistory
)
513 mdof_flags
|= MDOF_X
;
517 mdof_flags
|= MDOF_F
;
520 /* If we want IMD output, set appropriate MDOF flag */
523 mdof_flags
|= MDOF_IMD
;
526 mdoutf_write_to_trajectory_files(fplog
, cr
, outf
, mdof_flags
,
527 top_global
, step
, (double)step
,
528 &state
->s
, state_global
, observablesHistory
,
531 if (confout
!= nullptr && MASTER(cr
))
533 GMX_RELEASE_ASSERT(bX
, "The code below assumes that (with domain decomposition), x is collected to state_global in the call above.");
534 /* With domain decomposition the call above collected the state->s.x
535 * into state_global->x. Without DD we copy the local state pointer.
537 if (!DOMAINDECOMP(cr
))
539 state_global
= &state
->s
;
542 if (ir
->ePBC
!= epbcNONE
&& !ir
->bPeriodicMols
&& DOMAINDECOMP(cr
))
544 /* Make molecules whole only for confout writing */
545 do_pbc_mtop(fplog
, ir
->ePBC
, state
->s
.box
, top_global
,
546 as_rvec_array(state_global
->x
.data()));
549 write_sto_conf_mtop(confout
,
550 *top_global
->name
, top_global
,
551 as_rvec_array(state_global
->x
.data()), nullptr, ir
->ePBC
, state
->s
.box
);
555 //! \brief Do one minimization step
557 // \returns true when the step succeeded, false when a constraint error occurred
558 static bool do_em_step(t_commrec
*cr
, t_inputrec
*ir
, t_mdatoms
*md
,
560 em_state_t
*ems1
, real a
, const PaddedRVecVector
*force
,
562 gmx_constr_t constr
, gmx_localtop_t
*top
,
563 t_nrnb
*nrnb
, gmx_wallcycle_t wcycle
,
570 int nthreads gmx_unused
;
572 bool validStep
= true;
577 if (DOMAINDECOMP(cr
) && s1
->ddp_count
!= cr
->dd
->ddp_count
)
579 gmx_incons("state mismatch in do_em_step");
582 s2
->flags
= s1
->flags
;
584 if (s2
->natoms
!= s1
->natoms
)
586 state_change_natoms(s2
, s1
->natoms
);
587 /* We need to allocate one element extra, since we might use
588 * (unaligned) 4-wide SIMD loads to access rvec entries.
590 ems2
->f
.resize(s2
->natoms
+ 1);
592 if (DOMAINDECOMP(cr
) && s2
->cg_gl
.size() != s1
->cg_gl
.size())
594 s2
->cg_gl
.resize(s1
->cg_gl
.size());
597 copy_mat(s1
->box
, s2
->box
);
598 /* Copy free energy state */
599 s2
->lambda
= s1
->lambda
;
600 copy_mat(s1
->box
, s2
->box
);
605 // cppcheck-suppress unreadVariable
606 nthreads
= gmx_omp_nthreads_get(emntUpdate
);
607 #pragma omp parallel num_threads(nthreads)
609 const rvec
*x1
= as_rvec_array(s1
->x
.data());
610 rvec
*x2
= as_rvec_array(s2
->x
.data());
611 const rvec
*f
= as_rvec_array(force
->data());
614 #pragma omp for schedule(static) nowait
615 for (int i
= start
; i
< end
; i
++)
623 for (int m
= 0; m
< DIM
; m
++)
625 if (ir
->opts
.nFreeze
[gf
][m
])
631 x2
[i
][m
] = x1
[i
][m
] + a
*f
[i
][m
];
635 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
638 if (s2
->flags
& (1<<estCGP
))
640 /* Copy the CG p vector */
641 const rvec
*p1
= as_rvec_array(s1
->cg_p
.data());
642 rvec
*p2
= as_rvec_array(s2
->cg_p
.data());
643 #pragma omp for schedule(static) nowait
644 for (int i
= start
; i
< end
; i
++)
646 // Trivial OpenMP block that does not throw
647 copy_rvec(p1
[i
], p2
[i
]);
651 if (DOMAINDECOMP(cr
))
653 s2
->ddp_count
= s1
->ddp_count
;
655 /* OpenMP does not supported unsigned loop variables */
656 #pragma omp for schedule(static) nowait
657 for (int i
= 0; i
< static_cast<int>(s2
->cg_gl
.size()); i
++)
659 s2
->cg_gl
[i
] = s1
->cg_gl
[i
];
661 s2
->ddp_count_cg_gl
= s1
->ddp_count_cg_gl
;
667 wallcycle_start(wcycle
, ewcCONSTR
);
670 constrain(nullptr, TRUE
, TRUE
, constr
, &top
->idef
,
671 ir
, cr
, count
, 0, 1.0, md
,
672 as_rvec_array(s1
->x
.data()), as_rvec_array(s2
->x
.data()),
673 nullptr, bMolPBC
, s2
->box
,
674 s2
->lambda
[efptBONDED
], &dvdl_constr
,
675 nullptr, nullptr, nrnb
, econqCoord
);
676 wallcycle_stop(wcycle
, ewcCONSTR
);
678 // We should move this check to the different minimizers
679 if (!validStep
&& ir
->eI
!= eiSteep
)
681 gmx_fatal(FARGS
, "The coordinates could not be constrained. Minimizer '%s' can not handle constraint failures, use minimizer '%s' before using '%s'.",
682 EI(ir
->eI
), EI(eiSteep
), EI(ir
->eI
));
689 //! Prepare EM for using domain decomposition parallellization
690 static void em_dd_partition_system(FILE *fplog
, int step
, t_commrec
*cr
,
691 gmx_mtop_t
*top_global
, t_inputrec
*ir
,
692 em_state_t
*ems
, gmx_localtop_t
*top
,
693 t_mdatoms
*mdatoms
, t_forcerec
*fr
,
694 gmx_vsite_t
*vsite
, gmx_constr_t constr
,
695 t_nrnb
*nrnb
, gmx_wallcycle_t wcycle
)
697 /* Repartition the domain decomposition */
698 dd_partition_system(fplog
, step
, cr
, FALSE
, 1,
699 nullptr, top_global
, ir
,
701 mdatoms
, top
, fr
, vsite
, constr
,
702 nrnb
, wcycle
, FALSE
);
703 dd_store_state(cr
->dd
, &ems
->s
);
706 //! De one energy evaluation
707 static void evaluate_energy(FILE *fplog
, t_commrec
*cr
,
708 gmx_mtop_t
*top_global
,
709 em_state_t
*ems
, gmx_localtop_t
*top
,
710 t_inputrec
*inputrec
,
711 t_nrnb
*nrnb
, gmx_wallcycle_t wcycle
,
712 gmx_global_stat_t gstat
,
713 gmx_vsite_t
*vsite
, gmx_constr_t constr
,
715 t_graph
*graph
, t_mdatoms
*mdatoms
,
716 t_forcerec
*fr
, rvec mu_tot
,
717 gmx_enerdata_t
*enerd
, tensor vir
, tensor pres
,
718 gmx_int64_t count
, gmx_bool bFirst
)
722 tensor force_vir
, shake_vir
, ekin
;
723 real dvdl_constr
, prescorr
, enercorr
, dvdlcorr
;
726 /* Set the time to the initial time, the time does not change during EM */
727 t
= inputrec
->init_t
;
730 (DOMAINDECOMP(cr
) && ems
->s
.ddp_count
< cr
->dd
->ddp_count
))
732 /* This is the first state or an old state used before the last ns */
738 if (inputrec
->nstlist
> 0)
746 construct_vsites(vsite
, as_rvec_array(ems
->s
.x
.data()), 1, nullptr,
747 top
->idef
.iparams
, top
->idef
.il
,
748 fr
->ePBC
, fr
->bMolPBC
, cr
, ems
->s
.box
);
751 if (DOMAINDECOMP(cr
) && bNS
)
753 /* Repartition the domain decomposition */
754 em_dd_partition_system(fplog
, count
, cr
, top_global
, inputrec
,
755 ems
, top
, mdatoms
, fr
, vsite
, constr
,
759 /* Calc force & energy on new trial position */
760 /* do_force always puts the charge groups in the box and shifts again
761 * We do not unshift, so molecules are always whole in congrad.c
763 do_force(fplog
, cr
, inputrec
,
764 count
, nrnb
, wcycle
, top
, &top_global
->groups
,
765 ems
->s
.box
, &ems
->s
.x
, &ems
->s
.hist
,
766 &ems
->f
, force_vir
, mdatoms
, enerd
, fcd
,
767 ems
->s
.lambda
, graph
, fr
, vsite
, mu_tot
, t
, nullptr, TRUE
,
768 GMX_FORCE_STATECHANGED
| GMX_FORCE_ALLFORCES
|
769 GMX_FORCE_VIRIAL
| GMX_FORCE_ENERGY
|
770 (bNS
? GMX_FORCE_NS
: 0));
772 /* Clear the unused shake virial and pressure */
773 clear_mat(shake_vir
);
776 /* Communicate stuff when parallel */
777 if (PAR(cr
) && inputrec
->eI
!= eiNM
)
779 wallcycle_start(wcycle
, ewcMoveE
);
781 global_stat(gstat
, cr
, enerd
, force_vir
, shake_vir
, mu_tot
,
782 inputrec
, nullptr, nullptr, nullptr, 1, &terminate
,
788 wallcycle_stop(wcycle
, ewcMoveE
);
791 /* Calculate long range corrections to pressure and energy */
792 calc_dispcorr(inputrec
, fr
, ems
->s
.box
, ems
->s
.lambda
[efptVDW
],
793 pres
, force_vir
, &prescorr
, &enercorr
, &dvdlcorr
);
794 enerd
->term
[F_DISPCORR
] = enercorr
;
795 enerd
->term
[F_EPOT
] += enercorr
;
796 enerd
->term
[F_PRES
] += prescorr
;
797 enerd
->term
[F_DVDL
] += dvdlcorr
;
799 ems
->epot
= enerd
->term
[F_EPOT
];
803 /* Project out the constraint components of the force */
804 wallcycle_start(wcycle
, ewcCONSTR
);
806 rvec
*f_rvec
= as_rvec_array(ems
->f
.data());
807 constrain(nullptr, FALSE
, FALSE
, constr
, &top
->idef
,
808 inputrec
, cr
, count
, 0, 1.0, mdatoms
,
809 as_rvec_array(ems
->s
.x
.data()), f_rvec
, f_rvec
,
810 fr
->bMolPBC
, ems
->s
.box
,
811 ems
->s
.lambda
[efptBONDED
], &dvdl_constr
,
812 nullptr, &shake_vir
, nrnb
, econqForceDispl
);
813 enerd
->term
[F_DVDL_CONSTR
] += dvdl_constr
;
814 m_add(force_vir
, shake_vir
, vir
);
815 wallcycle_stop(wcycle
, ewcCONSTR
);
819 copy_mat(force_vir
, vir
);
823 enerd
->term
[F_PRES
] =
824 calc_pres(fr
->ePBC
, inputrec
->nwall
, ems
->s
.box
, ekin
, vir
, pres
);
826 sum_dhdl(enerd
, ems
->s
.lambda
, inputrec
->fepvals
);
828 if (EI_ENERGY_MINIMIZATION(inputrec
->eI
))
830 get_state_f_norm_max(cr
, &(inputrec
->opts
), mdatoms
, ems
);
834 //! Parallel utility summing energies and forces
835 static double reorder_partsum(t_commrec
*cr
, t_grpopts
*opts
, t_mdatoms
*mdatoms
,
836 gmx_mtop_t
*top_global
,
837 em_state_t
*s_min
, em_state_t
*s_b
)
840 int ncg
, *cg_gl
, *index
, c
, cg
, i
, a0
, a1
, a
, gf
, m
;
842 unsigned char *grpnrFREEZE
;
846 fprintf(debug
, "Doing reorder_partsum\n");
849 const rvec
*fm
= as_rvec_array(s_min
->f
.data());
850 const rvec
*fb
= as_rvec_array(s_b
->f
.data());
852 cgs_gl
= dd_charge_groups_global(cr
->dd
);
853 index
= cgs_gl
->index
;
855 /* Collect fm in a global vector fmg.
856 * This conflicts with the spirit of domain decomposition,
857 * but to fully optimize this a much more complicated algorithm is required.
860 snew(fmg
, top_global
->natoms
);
862 ncg
= s_min
->s
.cg_gl
.size();
863 cg_gl
= s_min
->s
.cg_gl
.data();
865 for (c
= 0; c
< ncg
; c
++)
870 for (a
= a0
; a
< a1
; a
++)
872 copy_rvec(fm
[i
], fmg
[a
]);
876 gmx_sum(top_global
->natoms
*3, fmg
[0], cr
);
878 /* Now we will determine the part of the sum for the cgs in state s_b */
879 ncg
= s_b
->s
.cg_gl
.size();
880 cg_gl
= s_b
->s
.cg_gl
.data();
884 grpnrFREEZE
= top_global
->groups
.grpnr
[egcFREEZE
];
885 for (c
= 0; c
< ncg
; c
++)
890 for (a
= a0
; a
< a1
; a
++)
892 if (mdatoms
->cFREEZE
&& grpnrFREEZE
)
896 for (m
= 0; m
< DIM
; m
++)
898 if (!opts
->nFreeze
[gf
][m
])
900 partsum
+= (fb
[i
][m
] - fmg
[a
][m
])*fb
[i
][m
];
912 //! Print some stuff, like beta, whatever that means.
913 static real
pr_beta(t_commrec
*cr
, t_grpopts
*opts
, t_mdatoms
*mdatoms
,
914 gmx_mtop_t
*top_global
,
915 em_state_t
*s_min
, em_state_t
*s_b
)
919 /* This is just the classical Polak-Ribiere calculation of beta;
920 * it looks a bit complicated since we take freeze groups into account,
921 * and might have to sum it in parallel runs.
924 if (!DOMAINDECOMP(cr
) ||
925 (s_min
->s
.ddp_count
== cr
->dd
->ddp_count
&&
926 s_b
->s
.ddp_count
== cr
->dd
->ddp_count
))
928 const rvec
*fm
= as_rvec_array(s_min
->f
.data());
929 const rvec
*fb
= as_rvec_array(s_b
->f
.data());
932 /* This part of code can be incorrect with DD,
933 * since the atom ordering in s_b and s_min might differ.
935 for (int i
= 0; i
< mdatoms
->homenr
; i
++)
937 if (mdatoms
->cFREEZE
)
939 gf
= mdatoms
->cFREEZE
[i
];
941 for (int m
= 0; m
< DIM
; m
++)
943 if (!opts
->nFreeze
[gf
][m
])
945 sum
+= (fb
[i
][m
] - fm
[i
][m
])*fb
[i
][m
];
952 /* We need to reorder cgs while summing */
953 sum
= reorder_partsum(cr
, opts
, mdatoms
, top_global
, s_min
, s_b
);
957 gmx_sumd(1, &sum
, cr
);
960 return sum
/gmx::square(s_min
->fnorm
);
966 /*! \brief Do conjugate gradients minimization
967 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
968 int nfile, const t_filenm fnm[],
969 const gmx_output_env_t *oenv, gmx_bool bVerbose,
971 gmx_vsite_t *vsite, gmx_constr_t constr,
973 gmx::IMDOutputProvider *outputProvider,
974 t_inputrec *inputrec,
975 gmx_mtop_t *top_global, t_fcdata *fcd,
976 t_state *state_global,
978 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
981 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
982 gmx_membed_t gmx_unused *membed,
983 real cpt_period, real max_hours,
986 gmx_walltime_accounting_t walltime_accounting)
988 double do_cg(FILE *fplog
, t_commrec
*cr
, const gmx::MDLogger gmx_unused
&mdlog
,
989 int nfile
, const t_filenm fnm
[],
990 const gmx_output_env_t gmx_unused
*oenv
, gmx_bool bVerbose
,
991 int gmx_unused nstglobalcomm
,
992 gmx_vsite_t
*vsite
, gmx_constr_t constr
,
993 int gmx_unused stepout
,
994 gmx::IMDOutputProvider
*outputProvider
,
995 t_inputrec
*inputrec
,
996 gmx_mtop_t
*top_global
, t_fcdata
*fcd
,
997 t_state
*state_global
,
998 ObservablesHistory
*observablesHistory
,
1000 t_nrnb
*nrnb
, gmx_wallcycle_t wcycle
,
1002 int gmx_unused repl_ex_nst
, int gmx_unused repl_ex_nex
, int gmx_unused repl_ex_seed
,
1003 gmx_membed_t gmx_unused
*membed
,
1004 real gmx_unused cpt_period
, real gmx_unused max_hours
,
1006 unsigned long gmx_unused Flags
,
1007 gmx_walltime_accounting_t walltime_accounting
)
1009 const char *CG
= "Polak-Ribiere Conjugate Gradients";
1011 gmx_localtop_t
*top
;
1012 gmx_enerdata_t
*enerd
;
1013 gmx_global_stat_t gstat
;
1015 double tmp
, minstep
;
1017 real a
, b
, c
, beta
= 0.0;
1021 gmx_bool converged
, foundlower
;
1023 gmx_bool do_log
= FALSE
, do_ene
= FALSE
, do_x
, do_f
;
1025 int number_steps
, neval
= 0, nstcg
= inputrec
->nstcgsteep
;
1027 int m
, step
, nminstep
;
1031 // Ensure the extra per-atom state array gets allocated
1032 state_global
->flags
|= (1<<estCGP
);
1034 /* Create 4 states on the stack and extract pointers that we will swap */
1035 em_state_t s0
{}, s1
{}, s2
{}, s3
{};
1036 em_state_t
*s_min
= &s0
;
1037 em_state_t
*s_a
= &s1
;
1038 em_state_t
*s_b
= &s2
;
1039 em_state_t
*s_c
= &s3
;
1041 /* Init em and store the local state in s_min */
1042 init_em(fplog
, CG
, cr
, outputProvider
, inputrec
,
1043 state_global
, top_global
, s_min
, &top
,
1044 nrnb
, mu_tot
, fr
, &enerd
, &graph
, mdatoms
, &gstat
,
1045 vsite
, constr
, nullptr,
1046 nfile
, fnm
, &outf
, &mdebin
, imdport
, Flags
, wcycle
);
1048 /* Print to log file */
1049 print_em_start(fplog
, cr
, walltime_accounting
, wcycle
, CG
);
1051 /* Max number of steps */
1052 number_steps
= inputrec
->nsteps
;
1056 sp_header(stderr
, CG
, inputrec
->em_tol
, number_steps
);
1060 sp_header(fplog
, CG
, inputrec
->em_tol
, number_steps
);
1063 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1064 /* do_force always puts the charge groups in the box and shifts again
1065 * We do not unshift, so molecules are always whole in congrad.c
1067 evaluate_energy(fplog
, cr
,
1068 top_global
, s_min
, top
,
1069 inputrec
, nrnb
, wcycle
, gstat
,
1070 vsite
, constr
, fcd
, graph
, mdatoms
, fr
,
1071 mu_tot
, enerd
, vir
, pres
, -1, TRUE
);
1076 /* Copy stuff to the energy bin for easy printing etc. */
1077 upd_mdebin(mdebin
, FALSE
, FALSE
, (double)step
,
1078 mdatoms
->tmass
, enerd
, &s_min
->s
, inputrec
->fepvals
, inputrec
->expandedvals
, s_min
->s
.box
,
1079 nullptr, nullptr, vir
, pres
, nullptr, mu_tot
, constr
);
1081 print_ebin_header(fplog
, step
, step
);
1082 print_ebin(mdoutf_get_fp_ene(outf
), TRUE
, FALSE
, FALSE
, fplog
, step
, step
, eprNORMAL
,
1083 mdebin
, fcd
, &(top_global
->groups
), &(inputrec
->opts
));
1087 /* Estimate/guess the initial stepsize */
1088 stepsize
= inputrec
->em_stepsize
/s_min
->fnorm
;
1092 double sqrtNumAtoms
= sqrt(static_cast<double>(state_global
->natoms
));
1093 fprintf(stderr
, " F-max = %12.5e on atom %d\n",
1094 s_min
->fmax
, s_min
->a_fmax
+1);
1095 fprintf(stderr
, " F-Norm = %12.5e\n",
1096 s_min
->fnorm
/sqrtNumAtoms
);
1097 fprintf(stderr
, "\n");
1098 /* and copy to the log file too... */
1099 fprintf(fplog
, " F-max = %12.5e on atom %d\n",
1100 s_min
->fmax
, s_min
->a_fmax
+1);
1101 fprintf(fplog
, " F-Norm = %12.5e\n",
1102 s_min
->fnorm
/sqrtNumAtoms
);
1103 fprintf(fplog
, "\n");
1105 /* Start the loop over CG steps.
1106 * Each successful step is counted, and we continue until
1107 * we either converge or reach the max number of steps.
1110 for (step
= 0; (number_steps
< 0 || step
<= number_steps
) && !converged
; step
++)
1113 /* start taking steps in a new direction
1114 * First time we enter the routine, beta=0, and the direction is
1115 * simply the negative gradient.
1118 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1119 rvec
*pm
= as_rvec_array(s_min
->s
.cg_p
.data());
1120 const rvec
*sfm
= as_rvec_array(s_min
->f
.data());
1123 for (int i
= 0; i
< mdatoms
->homenr
; i
++)
1125 if (mdatoms
->cFREEZE
)
1127 gf
= mdatoms
->cFREEZE
[i
];
1129 for (m
= 0; m
< DIM
; m
++)
1131 if (!inputrec
->opts
.nFreeze
[gf
][m
])
1133 pm
[i
][m
] = sfm
[i
][m
] + beta
*pm
[i
][m
];
1134 gpa
-= pm
[i
][m
]*sfm
[i
][m
];
1135 /* f is negative gradient, thus the sign */
1144 /* Sum the gradient along the line across CPUs */
1147 gmx_sumd(1, &gpa
, cr
);
1150 /* Calculate the norm of the search vector */
1151 get_f_norm_max(cr
, &(inputrec
->opts
), mdatoms
, pm
, &pnorm
, nullptr, nullptr);
1153 /* Just in case stepsize reaches zero due to numerical precision... */
1156 stepsize
= inputrec
->em_stepsize
/pnorm
;
1160 * Double check the value of the derivative in the search direction.
1161 * If it is positive it must be due to the old information in the
1162 * CG formula, so just remove that and start over with beta=0.
1163 * This corresponds to a steepest descent step.
1168 step
--; /* Don't count this step since we are restarting */
1169 continue; /* Go back to the beginning of the big for-loop */
1172 /* Calculate minimum allowed stepsize, before the average (norm)
1173 * relative change in coordinate is smaller than precision
1176 for (int i
= 0; i
< mdatoms
->homenr
; i
++)
1178 for (m
= 0; m
< DIM
; m
++)
1180 tmp
= fabs(s_min
->s
.x
[i
][m
]);
1189 /* Add up from all CPUs */
1192 gmx_sumd(1, &minstep
, cr
);
1195 minstep
= GMX_REAL_EPS
/sqrt(minstep
/(3*state_global
->natoms
));
1197 if (stepsize
< minstep
)
1203 /* Write coordinates if necessary */
1204 do_x
= do_per_step(step
, inputrec
->nstxout
);
1205 do_f
= do_per_step(step
, inputrec
->nstfout
);
1207 write_em_traj(fplog
, cr
, outf
, do_x
, do_f
, nullptr,
1208 top_global
, inputrec
, step
,
1209 s_min
, state_global
, observablesHistory
);
1211 /* Take a step downhill.
1212 * In theory, we should minimize the function along this direction.
1213 * That is quite possible, but it turns out to take 5-10 function evaluations
1214 * for each line. However, we dont really need to find the exact minimum -
1215 * it is much better to start a new CG step in a modified direction as soon
1216 * as we are close to it. This will save a lot of energy evaluations.
1218 * In practice, we just try to take a single step.
1219 * If it worked (i.e. lowered the energy), we increase the stepsize but
1220 * the continue straight to the next CG step without trying to find any minimum.
1221 * If it didn't work (higher energy), there must be a minimum somewhere between
1222 * the old position and the new one.
1224 * Due to the finite numerical accuracy, it turns out that it is a good idea
1225 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1226 * This leads to lower final energies in the tests I've done. / Erik
1228 s_a
->epot
= s_min
->epot
;
1230 c
= a
+ stepsize
; /* reference position along line is zero */
1232 if (DOMAINDECOMP(cr
) && s_min
->s
.ddp_count
< cr
->dd
->ddp_count
)
1234 em_dd_partition_system(fplog
, step
, cr
, top_global
, inputrec
,
1235 s_min
, top
, mdatoms
, fr
, vsite
, constr
,
1239 /* Take a trial step (new coords in s_c) */
1240 do_em_step(cr
, inputrec
, mdatoms
, fr
->bMolPBC
, s_min
, c
, &s_min
->s
.cg_p
, s_c
,
1241 constr
, top
, nrnb
, wcycle
, -1);
1244 /* Calculate energy for the trial step */
1245 evaluate_energy(fplog
, cr
,
1246 top_global
, s_c
, top
,
1247 inputrec
, nrnb
, wcycle
, gstat
,
1248 vsite
, constr
, fcd
, graph
, mdatoms
, fr
,
1249 mu_tot
, enerd
, vir
, pres
, -1, FALSE
);
1251 /* Calc derivative along line */
1252 const rvec
*pc
= as_rvec_array(s_c
->s
.cg_p
.data());
1253 const rvec
*sfc
= as_rvec_array(s_c
->f
.data());
1255 for (int i
= 0; i
< mdatoms
->homenr
; i
++)
1257 for (m
= 0; m
< DIM
; m
++)
1259 gpc
-= pc
[i
][m
]*sfc
[i
][m
]; /* f is negative gradient, thus the sign */
1262 /* Sum the gradient along the line across CPUs */
1265 gmx_sumd(1, &gpc
, cr
);
1268 /* This is the max amount of increase in energy we tolerate */
1269 tmp
= sqrt(GMX_REAL_EPS
)*fabs(s_a
->epot
);
1271 /* Accept the step if the energy is lower, or if it is not significantly higher
1272 * and the line derivative is still negative.
1274 if (s_c
->epot
< s_a
->epot
|| (gpc
< 0 && s_c
->epot
< (s_a
->epot
+ tmp
)))
1277 /* Great, we found a better energy. Increase step for next iteration
1278 * if we are still going down, decrease it otherwise
1282 stepsize
*= 1.618034; /* The golden section */
1286 stepsize
*= 0.618034; /* 1/golden section */
1291 /* New energy is the same or higher. We will have to do some work
1292 * to find a smaller value in the interval. Take smaller step next time!
1295 stepsize
*= 0.618034;
1301 /* OK, if we didn't find a lower value we will have to locate one now - there must
1302 * be one in the interval [a=0,c].
1303 * The same thing is valid here, though: Don't spend dozens of iterations to find
1304 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1305 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1307 * I also have a safeguard for potentially really pathological functions so we never
1308 * take more than 20 steps before we give up ...
1310 * If we already found a lower value we just skip this step and continue to the update.
1319 /* Select a new trial point.
1320 * If the derivatives at points a & c have different sign we interpolate to zero,
1321 * otherwise just do a bisection.
1323 if (gpa
< 0 && gpc
> 0)
1325 b
= a
+ gpa
*(a
-c
)/(gpc
-gpa
);
1332 /* safeguard if interpolation close to machine accuracy causes errors:
1333 * never go outside the interval
1335 if (b
<= a
|| b
>= c
)
1340 if (DOMAINDECOMP(cr
) && s_min
->s
.ddp_count
!= cr
->dd
->ddp_count
)
1342 /* Reload the old state */
1343 em_dd_partition_system(fplog
, -1, cr
, top_global
, inputrec
,
1344 s_min
, top
, mdatoms
, fr
, vsite
, constr
,
1348 /* Take a trial step to this new point - new coords in s_b */
1349 do_em_step(cr
, inputrec
, mdatoms
, fr
->bMolPBC
, s_min
, b
, &s_min
->s
.cg_p
, s_b
,
1350 constr
, top
, nrnb
, wcycle
, -1);
1353 /* Calculate energy for the trial step */
1354 evaluate_energy(fplog
, cr
,
1355 top_global
, s_b
, top
,
1356 inputrec
, nrnb
, wcycle
, gstat
,
1357 vsite
, constr
, fcd
, graph
, mdatoms
, fr
,
1358 mu_tot
, enerd
, vir
, pres
, -1, FALSE
);
1360 /* p does not change within a step, but since the domain decomposition
1361 * might change, we have to use cg_p of s_b here.
1363 const rvec
*pb
= as_rvec_array(s_b
->s
.cg_p
.data());
1364 const rvec
*sfb
= as_rvec_array(s_b
->f
.data());
1366 for (int i
= 0; i
< mdatoms
->homenr
; i
++)
1368 for (m
= 0; m
< DIM
; m
++)
1370 gpb
-= pb
[i
][m
]*sfb
[i
][m
]; /* f is negative gradient, thus the sign */
1373 /* Sum the gradient along the line across CPUs */
1376 gmx_sumd(1, &gpb
, cr
);
1381 fprintf(debug
, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1382 s_a
->epot
, s_b
->epot
, s_c
->epot
, gpb
);
1385 epot_repl
= s_b
->epot
;
1387 /* Keep one of the intervals based on the value of the derivative at the new point */
1390 /* Replace c endpoint with b */
1391 swap_em_state(&s_b
, &s_c
);
1397 /* Replace a endpoint with b */
1398 swap_em_state(&s_b
, &s_a
);
1404 * Stop search as soon as we find a value smaller than the endpoints.
1405 * Never run more than 20 steps, no matter what.
1409 while ((epot_repl
> s_a
->epot
|| epot_repl
> s_c
->epot
) &&
1412 if (fabs(epot_repl
- s_min
->epot
) < fabs(s_min
->epot
)*GMX_REAL_EPS
||
1415 /* OK. We couldn't find a significantly lower energy.
1416 * If beta==0 this was steepest descent, and then we give up.
1417 * If not, set beta=0 and restart with steepest descent before quitting.
1427 /* Reset memory before giving up */
1433 /* Select min energy state of A & C, put the best in B.
1435 if (s_c
->epot
< s_a
->epot
)
1439 fprintf(debug
, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1440 s_c
->epot
, s_a
->epot
);
1442 swap_em_state(&s_b
, &s_c
);
1449 fprintf(debug
, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1450 s_a
->epot
, s_c
->epot
);
1452 swap_em_state(&s_b
, &s_a
);
1461 fprintf(debug
, "CGE: Found a lower energy %f, moving C to B\n",
1464 swap_em_state(&s_b
, &s_c
);
1468 /* new search direction */
1469 /* beta = 0 means forget all memory and restart with steepest descents. */
1470 if (nstcg
&& ((step
% nstcg
) == 0))
1476 /* s_min->fnorm cannot be zero, because then we would have converged
1480 /* Polak-Ribiere update.
1481 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1483 beta
= pr_beta(cr
, &inputrec
->opts
, mdatoms
, top_global
, s_min
, s_b
);
1485 /* Limit beta to prevent oscillations */
1486 if (fabs(beta
) > 5.0)
1492 /* update positions */
1493 swap_em_state(&s_min
, &s_b
);
1496 /* Print it if necessary */
1501 double sqrtNumAtoms
= sqrt(static_cast<double>(state_global
->natoms
));
1502 fprintf(stderr
, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1503 step
, s_min
->epot
, s_min
->fnorm
/sqrtNumAtoms
,
1504 s_min
->fmax
, s_min
->a_fmax
+1);
1507 /* Store the new (lower) energies */
1508 upd_mdebin(mdebin
, FALSE
, FALSE
, (double)step
,
1509 mdatoms
->tmass
, enerd
, &s_min
->s
, inputrec
->fepvals
, inputrec
->expandedvals
, s_min
->s
.box
,
1510 nullptr, nullptr, vir
, pres
, nullptr, mu_tot
, constr
);
1512 do_log
= do_per_step(step
, inputrec
->nstlog
);
1513 do_ene
= do_per_step(step
, inputrec
->nstenergy
);
1515 /* Prepare IMD energy record, if bIMD is TRUE. */
1516 IMD_fill_energy_record(inputrec
->bIMD
, inputrec
->imd
, enerd
, step
, TRUE
);
1520 print_ebin_header(fplog
, step
, step
);
1522 print_ebin(mdoutf_get_fp_ene(outf
), do_ene
, FALSE
, FALSE
,
1523 do_log
? fplog
: nullptr, step
, step
, eprNORMAL
,
1524 mdebin
, fcd
, &(top_global
->groups
), &(inputrec
->opts
));
1527 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1528 if (do_IMD(inputrec
->bIMD
, step
, cr
, TRUE
, state_global
->box
, as_rvec_array(state_global
->x
.data()), inputrec
, 0, wcycle
) && MASTER(cr
))
1530 IMD_send_positions(inputrec
->imd
);
1533 /* Stop when the maximum force lies below tolerance.
1534 * If we have reached machine precision, converged is already set to true.
1536 converged
= converged
|| (s_min
->fmax
< inputrec
->em_tol
);
1538 } /* End of the loop */
1540 /* IMD cleanup, if bIMD is TRUE. */
1541 IMD_finalize(inputrec
->bIMD
, inputrec
->imd
);
1545 step
--; /* we never took that last step in this case */
1548 if (s_min
->fmax
> inputrec
->em_tol
)
1552 warn_step(stderr
, inputrec
->em_tol
, step
-1 == number_steps
, FALSE
);
1553 warn_step(fplog
, inputrec
->em_tol
, step
-1 == number_steps
, FALSE
);
1560 /* If we printed energy and/or logfile last step (which was the last step)
1561 * we don't have to do it again, but otherwise print the final values.
1565 /* Write final value to log since we didn't do anything the last step */
1566 print_ebin_header(fplog
, step
, step
);
1568 if (!do_ene
|| !do_log
)
1570 /* Write final energy file entries */
1571 print_ebin(mdoutf_get_fp_ene(outf
), !do_ene
, FALSE
, FALSE
,
1572 !do_log
? fplog
: nullptr, step
, step
, eprNORMAL
,
1573 mdebin
, fcd
, &(top_global
->groups
), &(inputrec
->opts
));
1577 /* Print some stuff... */
1580 fprintf(stderr
, "\nwriting lowest energy coordinates.\n");
1584 * For accurate normal mode calculation it is imperative that we
1585 * store the last conformation into the full precision binary trajectory.
1587 * However, we should only do it if we did NOT already write this step
1588 * above (which we did if do_x or do_f was true).
1590 do_x
= !do_per_step(step
, inputrec
->nstxout
);
1591 do_f
= (inputrec
->nstfout
> 0 && !do_per_step(step
, inputrec
->nstfout
));
1593 write_em_traj(fplog
, cr
, outf
, do_x
, do_f
, ftp2fn(efSTO
, nfile
, fnm
),
1594 top_global
, inputrec
, step
,
1595 s_min
, state_global
, observablesHistory
);
1600 double sqrtNumAtoms
= sqrt(static_cast<double>(state_global
->natoms
));
1601 print_converged(stderr
, CG
, inputrec
->em_tol
, step
, converged
, number_steps
,
1602 s_min
, sqrtNumAtoms
);
1603 print_converged(fplog
, CG
, inputrec
->em_tol
, step
, converged
, number_steps
,
1604 s_min
, sqrtNumAtoms
);
1606 fprintf(fplog
, "\nPerformed %d energy evaluations in total.\n", neval
);
1609 finish_em(cr
, outf
, walltime_accounting
, wcycle
);
1611 /* To print the actual number of steps we needed somewhere */
1612 walltime_accounting_set_nsteps_done(walltime_accounting
, step
);
1615 } /* That's all folks */
1618 /*! \brief Do L-BFGS conjugate gradients minimization
1619 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
1620 int nfile, const t_filenm fnm[],
1621 const gmx_output_env_t *oenv, gmx_bool bVerbose,
1623 gmx_vsite_t *vsite, gmx_constr_t constr,
1625 gmx::IMDOutputProvider *outputProvider,
1626 t_inputrec *inputrec,
1627 gmx_mtop_t *top_global, t_fcdata *fcd,
1628 t_state *state_global,
1630 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1633 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
1634 gmx_membed_t gmx_unused *membed,
1635 real cpt_period, real max_hours,
1637 unsigned long Flags,
1638 gmx_walltime_accounting_t walltime_accounting)
1640 double do_lbfgs(FILE *fplog
, t_commrec
*cr
, const gmx::MDLogger gmx_unused
&mdlog
,
1641 int nfile
, const t_filenm fnm
[],
1642 const gmx_output_env_t gmx_unused
*oenv
, gmx_bool bVerbose
,
1643 int gmx_unused nstglobalcomm
,
1644 gmx_vsite_t
*vsite
, gmx_constr_t constr
,
1645 int gmx_unused stepout
,
1646 gmx::IMDOutputProvider
*outputProvider
,
1647 t_inputrec
*inputrec
,
1648 gmx_mtop_t
*top_global
, t_fcdata
*fcd
,
1649 t_state
*state_global
,
1650 ObservablesHistory
*observablesHistory
,
1652 t_nrnb
*nrnb
, gmx_wallcycle_t wcycle
,
1654 int gmx_unused repl_ex_nst
, int gmx_unused repl_ex_nex
, int gmx_unused repl_ex_seed
,
1655 gmx_membed_t gmx_unused
*membed
,
1656 real gmx_unused cpt_period
, real gmx_unused max_hours
,
1658 unsigned long gmx_unused Flags
,
1659 gmx_walltime_accounting_t walltime_accounting
)
1661 static const char *LBFGS
= "Low-Memory BFGS Minimizer";
1663 gmx_localtop_t
*top
;
1664 gmx_enerdata_t
*enerd
;
1665 gmx_global_stat_t gstat
;
1667 int ncorr
, nmaxcorr
, point
, cp
, neval
, nminstep
;
1668 double stepsize
, step_taken
, gpa
, gpb
, gpc
, tmp
, minstep
;
1669 real
*rho
, *alpha
, *p
, *s
, **dx
, **dg
;
1670 real a
, b
, c
, maxdelta
, delta
;
1672 real dgdx
, dgdg
, sq
, yr
, beta
;
1676 gmx_bool do_log
, do_ene
, do_x
, do_f
, foundlower
, *frozen
;
1678 int start
, end
, number_steps
;
1680 int i
, k
, m
, n
, gf
, step
;
1685 gmx_fatal(FARGS
, "Cannot do parallel L-BFGS Minimization - yet.\n");
1688 if (nullptr != constr
)
1690 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).");
1693 n
= 3*state_global
->natoms
;
1694 nmaxcorr
= inputrec
->nbfgscorr
;
1699 snew(rho
, nmaxcorr
);
1700 snew(alpha
, nmaxcorr
);
1703 for (i
= 0; i
< nmaxcorr
; i
++)
1709 for (i
= 0; i
< nmaxcorr
; i
++)
1718 init_em(fplog
, LBFGS
, cr
, outputProvider
, inputrec
,
1719 state_global
, top_global
, &ems
, &top
,
1720 nrnb
, mu_tot
, fr
, &enerd
, &graph
, mdatoms
, &gstat
,
1721 vsite
, constr
, nullptr,
1722 nfile
, fnm
, &outf
, &mdebin
, imdport
, Flags
, wcycle
);
1725 end
= mdatoms
->homenr
;
1727 /* We need 4 working states */
1728 em_state_t s0
{}, s1
{}, s2
{}, s3
{};
1729 em_state_t
*sa
= &s0
;
1730 em_state_t
*sb
= &s1
;
1731 em_state_t
*sc
= &s2
;
1732 em_state_t
*last
= &s3
;
1733 /* Initialize by copying the state from ems (we could skip x and f here) */
1738 /* Print to log file */
1739 print_em_start(fplog
, cr
, walltime_accounting
, wcycle
, LBFGS
);
1741 do_log
= do_ene
= do_x
= do_f
= TRUE
;
1743 /* Max number of steps */
1744 number_steps
= inputrec
->nsteps
;
1746 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1748 for (i
= start
; i
< end
; i
++)
1750 if (mdatoms
->cFREEZE
)
1752 gf
= mdatoms
->cFREEZE
[i
];
1754 for (m
= 0; m
< DIM
; m
++)
1756 frozen
[3*i
+m
] = inputrec
->opts
.nFreeze
[gf
][m
];
1761 sp_header(stderr
, LBFGS
, inputrec
->em_tol
, number_steps
);
1765 sp_header(fplog
, LBFGS
, inputrec
->em_tol
, number_steps
);
1770 construct_vsites(vsite
, as_rvec_array(state_global
->x
.data()), 1, nullptr,
1771 top
->idef
.iparams
, top
->idef
.il
,
1772 fr
->ePBC
, fr
->bMolPBC
, cr
, state_global
->box
);
1775 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1776 /* do_force always puts the charge groups in the box and shifts again
1777 * We do not unshift, so molecules are always whole
1780 evaluate_energy(fplog
, cr
,
1781 top_global
, &ems
, top
,
1782 inputrec
, nrnb
, wcycle
, gstat
,
1783 vsite
, constr
, fcd
, graph
, mdatoms
, fr
,
1784 mu_tot
, enerd
, vir
, pres
, -1, TRUE
);
1789 /* Copy stuff to the energy bin for easy printing etc. */
1790 upd_mdebin(mdebin
, FALSE
, FALSE
, (double)step
,
1791 mdatoms
->tmass
, enerd
, state_global
, inputrec
->fepvals
, inputrec
->expandedvals
, state_global
->box
,
1792 nullptr, nullptr, vir
, pres
, nullptr, mu_tot
, constr
);
1794 print_ebin_header(fplog
, step
, step
);
1795 print_ebin(mdoutf_get_fp_ene(outf
), TRUE
, FALSE
, FALSE
, fplog
, step
, step
, eprNORMAL
,
1796 mdebin
, fcd
, &(top_global
->groups
), &(inputrec
->opts
));
1800 /* Set the initial step.
1801 * since it will be multiplied by the non-normalized search direction
1802 * vector (force vector the first time), we scale it by the
1803 * norm of the force.
1808 double sqrtNumAtoms
= sqrt(static_cast<double>(state_global
->natoms
));
1809 fprintf(stderr
, "Using %d BFGS correction steps.\n\n", nmaxcorr
);
1810 fprintf(stderr
, " F-max = %12.5e on atom %d\n", ems
.fmax
, ems
.a_fmax
+ 1);
1811 fprintf(stderr
, " F-Norm = %12.5e\n", ems
.fnorm
/sqrtNumAtoms
);
1812 fprintf(stderr
, "\n");
1813 /* and copy to the log file too... */
1814 fprintf(fplog
, "Using %d BFGS correction steps.\n\n", nmaxcorr
);
1815 fprintf(fplog
, " F-max = %12.5e on atom %d\n", ems
.fmax
, ems
.a_fmax
+ 1);
1816 fprintf(fplog
, " F-Norm = %12.5e\n", ems
.fnorm
/sqrtNumAtoms
);
1817 fprintf(fplog
, "\n");
1820 // Point is an index to the memory of search directions, where 0 is the first one.
1823 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1824 real
*fInit
= static_cast<real
*>(as_rvec_array(ems
.f
.data())[0]);
1825 for (i
= 0; i
< n
; i
++)
1829 dx
[point
][i
] = fInit
[i
]; /* Initial search direction */
1837 // Stepsize will be modified during the search, and actually it is not critical
1838 // (the main efficiency in the algorithm comes from changing directions), but
1839 // we still need an initial value, so estimate it as the inverse of the norm
1840 // so we take small steps where the potential fluctuates a lot.
1841 stepsize
= 1.0/ems
.fnorm
;
1843 /* Start the loop over BFGS steps.
1844 * Each successful step is counted, and we continue until
1845 * we either converge or reach the max number of steps.
1850 /* Set the gradient from the force */
1852 for (step
= 0; (number_steps
< 0 || step
<= number_steps
) && !converged
; step
++)
1855 /* Write coordinates if necessary */
1856 do_x
= do_per_step(step
, inputrec
->nstxout
);
1857 do_f
= do_per_step(step
, inputrec
->nstfout
);
1862 mdof_flags
|= MDOF_X
;
1867 mdof_flags
|= MDOF_F
;
1872 mdof_flags
|= MDOF_IMD
;
1875 mdoutf_write_to_trajectory_files(fplog
, cr
, outf
, mdof_flags
,
1876 top_global
, step
, (real
)step
, &ems
.s
, state_global
, observablesHistory
, &ems
.f
);
1878 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1880 /* make s a pointer to current search direction - point=0 first time we get here */
1883 real
*xx
= static_cast<real
*>(as_rvec_array(ems
.s
.x
.data())[0]);
1884 real
*ff
= static_cast<real
*>(as_rvec_array(ems
.f
.data())[0]);
1886 // calculate line gradient in position A
1887 for (gpa
= 0, i
= 0; i
< n
; i
++)
1892 /* Calculate minimum allowed stepsize along the line, before the average (norm)
1893 * relative change in coordinate is smaller than precision
1895 for (minstep
= 0, i
= 0; i
< n
; i
++)
1905 minstep
= GMX_REAL_EPS
/sqrt(minstep
/n
);
1907 if (stepsize
< minstep
)
1913 // Before taking any steps along the line, store the old position
1915 real
*lastx
= static_cast<real
*>(as_rvec_array(last
->s
.x
.data())[0]);
1916 real
*lastf
= static_cast<real
*>(as_rvec_array(last
->f
.data())[0]);
1921 /* Take a step downhill.
1922 * In theory, we should find the actual minimum of the function in this
1923 * direction, somewhere along the line.
1924 * That is quite possible, but it turns out to take 5-10 function evaluations
1925 * for each line. However, we dont really need to find the exact minimum -
1926 * it is much better to start a new BFGS step in a modified direction as soon
1927 * as we are close to it. This will save a lot of energy evaluations.
1929 * In practice, we just try to take a single step.
1930 * If it worked (i.e. lowered the energy), we increase the stepsize but
1931 * continue straight to the next BFGS step without trying to find any minimum,
1932 * i.e. we change the search direction too. If the line was smooth, it is
1933 * likely we are in a smooth region, and then it makes sense to take longer
1934 * steps in the modified search direction too.
1936 * If it didn't work (higher energy), there must be a minimum somewhere between
1937 * the old position and the new one. Then we need to start by finding a lower
1938 * value before we change search direction. Since the energy was apparently
1939 * quite rough, we need to decrease the step size.
1941 * Due to the finite numerical accuracy, it turns out that it is a good idea
1942 * to accept a SMALL increase in energy, if the derivative is still downhill.
1943 * This leads to lower final energies in the tests I've done. / Erik
1946 // State "A" is the first position along the line.
1947 // reference position along line is initially zero
1950 // Check stepsize first. We do not allow displacements
1951 // larger than emstep.
1955 // Pick a new position C by adding stepsize to A.
1958 // Calculate what the largest change in any individual coordinate
1959 // would be (translation along line * gradient along line)
1961 for (i
= 0; i
< n
; i
++)
1964 if (delta
> maxdelta
)
1969 // If any displacement is larger than the stepsize limit, reduce the step
1970 if (maxdelta
> inputrec
->em_stepsize
)
1975 while (maxdelta
> inputrec
->em_stepsize
);
1977 // Take a trial step and move the coordinate array xc[] to position C
1978 real
*xc
= static_cast<real
*>(as_rvec_array(sc
->s
.x
.data())[0]);
1979 for (i
= 0; i
< n
; i
++)
1981 xc
[i
] = lastx
[i
] + c
*s
[i
];
1985 // Calculate energy for the trial step in position C
1986 evaluate_energy(fplog
, cr
,
1987 top_global
, sc
, top
,
1988 inputrec
, nrnb
, wcycle
, gstat
,
1989 vsite
, constr
, fcd
, graph
, mdatoms
, fr
,
1990 mu_tot
, enerd
, vir
, pres
, step
, FALSE
);
1992 // Calc line gradient in position C
1993 real
*fc
= static_cast<real
*>(as_rvec_array(sc
->f
.data())[0]);
1994 for (gpc
= 0, i
= 0; i
< n
; i
++)
1996 gpc
-= s
[i
]*fc
[i
]; /* f is negative gradient, thus the sign */
1998 /* Sum the gradient along the line across CPUs */
2001 gmx_sumd(1, &gpc
, cr
);
2004 // This is the max amount of increase in energy we tolerate.
2005 // By allowing VERY small changes (close to numerical precision) we
2006 // frequently find even better (lower) final energies.
2007 tmp
= sqrt(GMX_REAL_EPS
)*fabs(sa
->epot
);
2009 // Accept the step if the energy is lower in the new position C (compared to A),
2010 // or if it is not significantly higher and the line derivative is still negative.
2011 if (sc
->epot
< sa
->epot
|| (gpc
< 0 && sc
->epot
< (sa
->epot
+ tmp
)))
2013 // Great, we found a better energy. We no longer try to alter the
2014 // stepsize, but simply accept this new better position. The we select a new
2015 // search direction instead, which will be much more efficient than continuing
2016 // to take smaller steps along a line. Set fnorm based on the new C position,
2017 // which will be used to update the stepsize to 1/fnorm further down.
2022 // If we got here, the energy is NOT lower in point C, i.e. it will be the same
2023 // or higher than in point A. In this case it is pointless to move to point C,
2024 // so we will have to do more iterations along the same line to find a smaller
2025 // value in the interval [A=0.0,C].
2026 // Here, A is still 0.0, but that will change when we do a search in the interval
2027 // [0.0,C] below. That search we will do by interpolation or bisection rather
2028 // than with the stepsize, so no need to modify it. For the next search direction
2029 // it will be reset to 1/fnorm anyway.
2035 // OK, if we didn't find a lower value we will have to locate one now - there must
2036 // be one in the interval [a,c].
2037 // The same thing is valid here, though: Don't spend dozens of iterations to find
2038 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2039 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2040 // I also have a safeguard for potentially really pathological functions so we never
2041 // take more than 20 steps before we give up.
2042 // If we already found a lower value we just skip this step and continue to the update.
2047 // Select a new trial point B in the interval [A,C].
2048 // If the derivatives at points a & c have different sign we interpolate to zero,
2049 // otherwise just do a bisection since there might be multiple minima/maxima
2050 // inside the interval.
2051 if (gpa
< 0 && gpc
> 0)
2053 b
= a
+ gpa
*(a
-c
)/(gpc
-gpa
);
2060 /* safeguard if interpolation close to machine accuracy causes errors:
2061 * never go outside the interval
2063 if (b
<= a
|| b
>= c
)
2068 // Take a trial step to point B
2069 real
*xb
= static_cast<real
*>(as_rvec_array(sb
->s
.x
.data())[0]);
2070 for (i
= 0; i
< n
; i
++)
2072 xb
[i
] = lastx
[i
] + b
*s
[i
];
2076 // Calculate energy for the trial step in point B
2077 evaluate_energy(fplog
, cr
,
2078 top_global
, sb
, top
,
2079 inputrec
, nrnb
, wcycle
, gstat
,
2080 vsite
, constr
, fcd
, graph
, mdatoms
, fr
,
2081 mu_tot
, enerd
, vir
, pres
, step
, FALSE
);
2084 // Calculate gradient in point B
2085 real
*fb
= static_cast<real
*>(as_rvec_array(sb
->f
.data())[0]);
2086 for (gpb
= 0, i
= 0; i
< n
; i
++)
2088 gpb
-= s
[i
]*fb
[i
]; /* f is negative gradient, thus the sign */
2091 /* Sum the gradient along the line across CPUs */
2094 gmx_sumd(1, &gpb
, cr
);
2097 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2098 // at the new point B, and rename the endpoints of this new interval A and C.
2101 /* Replace c endpoint with b */
2103 /* swap states b and c */
2104 swap_em_state(&sb
, &sc
);
2108 /* Replace a endpoint with b */
2110 /* swap states a and b */
2111 swap_em_state(&sa
, &sb
);
2115 * Stop search as soon as we find a value smaller than the endpoints,
2116 * or if the tolerance is below machine precision.
2117 * Never run more than 20 steps, no matter what.
2121 while ((sb
->epot
> sa
->epot
|| sb
->epot
> sc
->epot
) && (nminstep
< 20));
2123 if (fabs(sb
->epot
- Epot0
) < GMX_REAL_EPS
|| nminstep
>= 20)
2125 /* OK. We couldn't find a significantly lower energy.
2126 * If ncorr==0 this was steepest descent, and then we give up.
2127 * If not, reset memory to restart as steepest descent before quitting.
2139 /* Search in gradient direction */
2140 for (i
= 0; i
< n
; i
++)
2142 dx
[point
][i
] = ff
[i
];
2144 /* Reset stepsize */
2145 stepsize
= 1.0/fnorm
;
2150 /* Select min energy state of A & C, put the best in xx/ff/Epot
2152 if (sc
->epot
< sa
->epot
)
2174 /* Update the memory information, and calculate a new
2175 * approximation of the inverse hessian
2178 /* Have new data in Epot, xx, ff */
2179 if (ncorr
< nmaxcorr
)
2184 for (i
= 0; i
< n
; i
++)
2186 dg
[point
][i
] = lastf
[i
]-ff
[i
];
2187 dx
[point
][i
] *= step_taken
;
2192 for (i
= 0; i
< n
; i
++)
2194 dgdg
+= dg
[point
][i
]*dg
[point
][i
];
2195 dgdx
+= dg
[point
][i
]*dx
[point
][i
];
2200 rho
[point
] = 1.0/dgdx
;
2203 if (point
>= nmaxcorr
)
2209 for (i
= 0; i
< n
; i
++)
2216 /* Recursive update. First go back over the memory points */
2217 for (k
= 0; k
< ncorr
; k
++)
2226 for (i
= 0; i
< n
; i
++)
2228 sq
+= dx
[cp
][i
]*p
[i
];
2231 alpha
[cp
] = rho
[cp
]*sq
;
2233 for (i
= 0; i
< n
; i
++)
2235 p
[i
] -= alpha
[cp
]*dg
[cp
][i
];
2239 for (i
= 0; i
< n
; i
++)
2244 /* And then go forward again */
2245 for (k
= 0; k
< ncorr
; k
++)
2248 for (i
= 0; i
< n
; i
++)
2250 yr
+= p
[i
]*dg
[cp
][i
];
2254 beta
= alpha
[cp
]-beta
;
2256 for (i
= 0; i
< n
; i
++)
2258 p
[i
] += beta
*dx
[cp
][i
];
2268 for (i
= 0; i
< n
; i
++)
2272 dx
[point
][i
] = p
[i
];
2280 /* Print it if necessary */
2285 double sqrtNumAtoms
= sqrt(static_cast<double>(state_global
->natoms
));
2286 fprintf(stderr
, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2287 step
, ems
.epot
, ems
.fnorm
/sqrtNumAtoms
, ems
.fmax
, ems
.a_fmax
+ 1);
2290 /* Store the new (lower) energies */
2291 upd_mdebin(mdebin
, FALSE
, FALSE
, (double)step
,
2292 mdatoms
->tmass
, enerd
, state_global
, inputrec
->fepvals
, inputrec
->expandedvals
, state_global
->box
,
2293 nullptr, nullptr, vir
, pres
, nullptr, mu_tot
, constr
);
2294 do_log
= do_per_step(step
, inputrec
->nstlog
);
2295 do_ene
= do_per_step(step
, inputrec
->nstenergy
);
2298 print_ebin_header(fplog
, step
, step
);
2300 print_ebin(mdoutf_get_fp_ene(outf
), do_ene
, FALSE
, FALSE
,
2301 do_log
? fplog
: nullptr, step
, step
, eprNORMAL
,
2302 mdebin
, fcd
, &(top_global
->groups
), &(inputrec
->opts
));
2305 /* Send x and E to IMD client, if bIMD is TRUE. */
2306 if (do_IMD(inputrec
->bIMD
, step
, cr
, TRUE
, state_global
->box
, as_rvec_array(state_global
->x
.data()), inputrec
, 0, wcycle
) && MASTER(cr
))
2308 IMD_send_positions(inputrec
->imd
);
2311 // Reset stepsize in we are doing more iterations
2312 stepsize
= 1.0/ems
.fnorm
;
2314 /* Stop when the maximum force lies below tolerance.
2315 * If we have reached machine precision, converged is already set to true.
2317 converged
= converged
|| (ems
.fmax
< inputrec
->em_tol
);
2319 } /* End of the loop */
2321 /* IMD cleanup, if bIMD is TRUE. */
2322 IMD_finalize(inputrec
->bIMD
, inputrec
->imd
);
2326 step
--; /* we never took that last step in this case */
2329 if (ems
.fmax
> inputrec
->em_tol
)
2333 warn_step(stderr
, inputrec
->em_tol
, step
-1 == number_steps
, FALSE
);
2334 warn_step(fplog
, inputrec
->em_tol
, step
-1 == number_steps
, FALSE
);
2339 /* If we printed energy and/or logfile last step (which was the last step)
2340 * we don't have to do it again, but otherwise print the final values.
2342 if (!do_log
) /* Write final value to log since we didn't do anythin last step */
2344 print_ebin_header(fplog
, step
, step
);
2346 if (!do_ene
|| !do_log
) /* Write final energy file entries */
2348 print_ebin(mdoutf_get_fp_ene(outf
), !do_ene
, FALSE
, FALSE
,
2349 !do_log
? fplog
: nullptr, step
, step
, eprNORMAL
,
2350 mdebin
, fcd
, &(top_global
->groups
), &(inputrec
->opts
));
2353 /* Print some stuff... */
2356 fprintf(stderr
, "\nwriting lowest energy coordinates.\n");
2360 * For accurate normal mode calculation it is imperative that we
2361 * store the last conformation into the full precision binary trajectory.
2363 * However, we should only do it if we did NOT already write this step
2364 * above (which we did if do_x or do_f was true).
2366 do_x
= !do_per_step(step
, inputrec
->nstxout
);
2367 do_f
= !do_per_step(step
, inputrec
->nstfout
);
2368 write_em_traj(fplog
, cr
, outf
, do_x
, do_f
, ftp2fn(efSTO
, nfile
, fnm
),
2369 top_global
, inputrec
, step
,
2370 &ems
, state_global
, observablesHistory
);
2374 double sqrtNumAtoms
= sqrt(static_cast<double>(state_global
->natoms
));
2375 print_converged(stderr
, LBFGS
, inputrec
->em_tol
, step
, converged
,
2376 number_steps
, &ems
, sqrtNumAtoms
);
2377 print_converged(fplog
, LBFGS
, inputrec
->em_tol
, step
, converged
,
2378 number_steps
, &ems
, sqrtNumAtoms
);
2380 fprintf(fplog
, "\nPerformed %d energy evaluations in total.\n", neval
);
2383 finish_em(cr
, outf
, walltime_accounting
, wcycle
);
2385 /* To print the actual number of steps we needed somewhere */
2386 walltime_accounting_set_nsteps_done(walltime_accounting
, step
);
2389 } /* That's all folks */
2391 /*! \brief Do steepest descents minimization
2392 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
2393 int nfile, const t_filenm fnm[],
2394 const gmx_output_env_t *oenv, gmx_bool bVerbose,
2396 gmx_vsite_t *vsite, gmx_constr_t constr,
2398 gmx::IMDOutputProvider *outputProvider,
2399 t_inputrec *inputrec,
2400 gmx_mtop_t *top_global, t_fcdata *fcd,
2401 t_state *state_global,
2403 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2406 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
2407 real cpt_period, real max_hours,
2409 unsigned long Flags,
2410 gmx_walltime_accounting_t walltime_accounting)
2412 double do_steep(FILE *fplog
, t_commrec
*cr
, const gmx::MDLogger gmx_unused
&mdlog
,
2413 int nfile
, const t_filenm fnm
[],
2414 const gmx_output_env_t gmx_unused
*oenv
, gmx_bool bVerbose
,
2415 int gmx_unused nstglobalcomm
,
2416 gmx_vsite_t
*vsite
, gmx_constr_t constr
,
2417 int gmx_unused stepout
,
2418 gmx::IMDOutputProvider
*outputProvider
,
2419 t_inputrec
*inputrec
,
2420 gmx_mtop_t
*top_global
, t_fcdata
*fcd
,
2421 t_state
*state_global
,
2422 ObservablesHistory
*observablesHistory
,
2424 t_nrnb
*nrnb
, gmx_wallcycle_t wcycle
,
2426 int gmx_unused repl_ex_nst
, int gmx_unused repl_ex_nex
, int gmx_unused repl_ex_seed
,
2427 gmx_membed_t gmx_unused
*membed
,
2428 real gmx_unused cpt_period
, real gmx_unused max_hours
,
2430 unsigned long gmx_unused Flags
,
2431 gmx_walltime_accounting_t walltime_accounting
)
2433 const char *SD
= "Steepest Descents";
2434 gmx_localtop_t
*top
;
2435 gmx_enerdata_t
*enerd
;
2436 gmx_global_stat_t gstat
;
2442 gmx_bool bDone
, bAbort
, do_x
, do_f
;
2447 int steps_accepted
= 0;
2449 /* Create 2 states on the stack and extract pointers that we will swap */
2450 em_state_t s0
{}, s1
{};
2451 em_state_t
*s_min
= &s0
;
2452 em_state_t
*s_try
= &s1
;
2454 /* Init em and store the local state in s_try */
2455 init_em(fplog
, SD
, cr
, outputProvider
, inputrec
,
2456 state_global
, top_global
, s_try
, &top
,
2457 nrnb
, mu_tot
, fr
, &enerd
, &graph
, mdatoms
, &gstat
,
2458 vsite
, constr
, nullptr,
2459 nfile
, fnm
, &outf
, &mdebin
, imdport
, Flags
, wcycle
);
2461 /* Print to log file */
2462 print_em_start(fplog
, cr
, walltime_accounting
, wcycle
, SD
);
2464 /* Set variables for stepsize (in nm). This is the largest
2465 * step that we are going to make in any direction.
2467 ustep
= inputrec
->em_stepsize
;
2470 /* Max number of steps */
2471 nsteps
= inputrec
->nsteps
;
2475 /* Print to the screen */
2476 sp_header(stderr
, SD
, inputrec
->em_tol
, nsteps
);
2480 sp_header(fplog
, SD
, inputrec
->em_tol
, nsteps
);
2483 /**** HERE STARTS THE LOOP ****
2484 * count is the counter for the number of steps
2485 * bDone will be TRUE when the minimization has converged
2486 * bAbort will be TRUE when nsteps steps have been performed or when
2487 * the stepsize becomes smaller than is reasonable for machine precision
2492 while (!bDone
&& !bAbort
)
2494 bAbort
= (nsteps
>= 0) && (count
== nsteps
);
2496 /* set new coordinates, except for first step */
2497 bool validStep
= true;
2501 do_em_step(cr
, inputrec
, mdatoms
, fr
->bMolPBC
,
2502 s_min
, stepsize
, &s_min
->f
, s_try
,
2503 constr
, top
, nrnb
, wcycle
, count
);
2508 evaluate_energy(fplog
, cr
,
2509 top_global
, s_try
, top
,
2510 inputrec
, nrnb
, wcycle
, gstat
,
2511 vsite
, constr
, fcd
, graph
, mdatoms
, fr
,
2512 mu_tot
, enerd
, vir
, pres
, count
, count
== 0);
2516 // Signal constraint error during stepping with energy=inf
2517 s_try
->epot
= std::numeric_limits
<real
>::infinity();
2522 print_ebin_header(fplog
, count
, count
);
2527 s_min
->epot
= s_try
->epot
;
2530 /* Print it if necessary */
2535 fprintf(stderr
, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2536 count
, ustep
, s_try
->epot
, s_try
->fmax
, s_try
->a_fmax
+1,
2537 ( (count
== 0) || (s_try
->epot
< s_min
->epot
) ) ? '\n' : '\r');
2541 if ( (count
== 0) || (s_try
->epot
< s_min
->epot
) )
2543 /* Store the new (lower) energies */
2544 upd_mdebin(mdebin
, FALSE
, FALSE
, (double)count
,
2545 mdatoms
->tmass
, enerd
, &s_try
->s
, inputrec
->fepvals
, inputrec
->expandedvals
,
2546 s_try
->s
.box
, nullptr, nullptr, vir
, pres
, nullptr, mu_tot
, constr
);
2548 /* Prepare IMD energy record, if bIMD is TRUE. */
2549 IMD_fill_energy_record(inputrec
->bIMD
, inputrec
->imd
, enerd
, count
, TRUE
);
2551 print_ebin(mdoutf_get_fp_ene(outf
), TRUE
,
2552 do_per_step(steps_accepted
, inputrec
->nstdisreout
),
2553 do_per_step(steps_accepted
, inputrec
->nstorireout
),
2554 fplog
, count
, count
, eprNORMAL
,
2555 mdebin
, fcd
, &(top_global
->groups
), &(inputrec
->opts
));
2560 /* Now if the new energy is smaller than the previous...
2561 * or if this is the first step!
2562 * or if we did random steps!
2565 if ( (count
== 0) || (s_try
->epot
< s_min
->epot
) )
2569 /* Test whether the convergence criterion is met... */
2570 bDone
= (s_try
->fmax
< inputrec
->em_tol
);
2572 /* Copy the arrays for force, positions and energy */
2573 /* The 'Min' array always holds the coords and forces of the minimal
2575 swap_em_state(&s_min
, &s_try
);
2581 /* Write to trn, if necessary */
2582 do_x
= do_per_step(steps_accepted
, inputrec
->nstxout
);
2583 do_f
= do_per_step(steps_accepted
, inputrec
->nstfout
);
2584 write_em_traj(fplog
, cr
, outf
, do_x
, do_f
, nullptr,
2585 top_global
, inputrec
, count
,
2586 s_min
, state_global
, observablesHistory
);
2590 /* If energy is not smaller make the step smaller... */
2593 if (DOMAINDECOMP(cr
) && s_min
->s
.ddp_count
!= cr
->dd
->ddp_count
)
2595 /* Reload the old state */
2596 em_dd_partition_system(fplog
, count
, cr
, top_global
, inputrec
,
2597 s_min
, top
, mdatoms
, fr
, vsite
, constr
,
2602 /* Determine new step */
2603 stepsize
= ustep
/s_min
->fmax
;
2605 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2607 if (count
== nsteps
|| ustep
< 1e-12)
2609 if (count
== nsteps
|| ustep
< 1e-6)
2614 warn_step(stderr
, inputrec
->em_tol
, count
== nsteps
, constr
!= nullptr);
2615 warn_step(fplog
, inputrec
->em_tol
, count
== nsteps
, constr
!= nullptr);
2620 /* Send IMD energies and positions, if bIMD is TRUE. */
2621 if (do_IMD(inputrec
->bIMD
, count
, cr
, TRUE
, state_global
->box
, as_rvec_array(state_global
->x
.data()), inputrec
, 0, wcycle
) && MASTER(cr
))
2623 IMD_send_positions(inputrec
->imd
);
2627 } /* End of the loop */
2629 /* IMD cleanup, if bIMD is TRUE. */
2630 IMD_finalize(inputrec
->bIMD
, inputrec
->imd
);
2632 /* Print some data... */
2635 fprintf(stderr
, "\nwriting lowest energy coordinates.\n");
2637 write_em_traj(fplog
, cr
, outf
, TRUE
, inputrec
->nstfout
, ftp2fn(efSTO
, nfile
, fnm
),
2638 top_global
, inputrec
, count
,
2639 s_min
, state_global
, observablesHistory
);
2643 double sqrtNumAtoms
= sqrt(static_cast<double>(state_global
->natoms
));
2645 print_converged(stderr
, SD
, inputrec
->em_tol
, count
, bDone
, nsteps
,
2646 s_min
, sqrtNumAtoms
);
2647 print_converged(fplog
, SD
, inputrec
->em_tol
, count
, bDone
, nsteps
,
2648 s_min
, sqrtNumAtoms
);
2651 finish_em(cr
, outf
, walltime_accounting
, wcycle
);
2653 /* To print the actual number of steps we needed somewhere */
2654 inputrec
->nsteps
= count
;
2656 walltime_accounting_set_nsteps_done(walltime_accounting
, count
);
2659 } /* That's all folks */
2661 /*! \brief Do normal modes analysis
2662 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
2663 int nfile, const t_filenm fnm[],
2664 const gmx_output_env_t *oenv, gmx_bool bVerbose,
2666 gmx_vsite_t *vsite, gmx_constr_t constr,
2668 gmx::IMDOutputProvider *outputProvider,
2669 t_inputrec *inputrec,
2670 gmx_mtop_t *top_global, t_fcdata *fcd,
2671 t_state *state_global,
2673 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2676 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
2677 real cpt_period, real max_hours,
2679 unsigned long Flags,
2680 gmx_walltime_accounting_t walltime_accounting)
2682 double do_nm(FILE *fplog
, t_commrec
*cr
, const gmx::MDLogger
&mdlog
,
2683 int nfile
, const t_filenm fnm
[],
2684 const gmx_output_env_t gmx_unused
*oenv
, gmx_bool bVerbose
,
2685 int gmx_unused nstglobalcomm
,
2686 gmx_vsite_t
*vsite
, gmx_constr_t constr
,
2687 int gmx_unused stepout
,
2688 gmx::IMDOutputProvider
*outputProvider
,
2689 t_inputrec
*inputrec
,
2690 gmx_mtop_t
*top_global
, t_fcdata
*fcd
,
2691 t_state
*state_global
,
2692 ObservablesHistory gmx_unused
*observablesHistory
,
2694 t_nrnb
*nrnb
, gmx_wallcycle_t wcycle
,
2696 int gmx_unused repl_ex_nst
, int gmx_unused repl_ex_nex
, int gmx_unused repl_ex_seed
,
2697 gmx_membed_t gmx_unused
*membed
,
2698 real gmx_unused cpt_period
, real gmx_unused max_hours
,
2700 unsigned long gmx_unused Flags
,
2701 gmx_walltime_accounting_t walltime_accounting
)
2703 const char *NM
= "Normal Mode Analysis";
2706 gmx_localtop_t
*top
;
2707 gmx_enerdata_t
*enerd
;
2708 gmx_global_stat_t gstat
;
2713 gmx_bool bSparse
; /* use sparse matrix storage format */
2715 gmx_sparsematrix_t
* sparse_matrix
= nullptr;
2716 real
* full_matrix
= nullptr;
2718 /* added with respect to mdrun */
2720 real der_range
= 10.0*sqrt(GMX_REAL_EPS
);
2722 bool bIsMaster
= MASTER(cr
);
2724 if (constr
!= nullptr)
2726 gmx_fatal(FARGS
, "Constraints present with Normal Mode Analysis, this combination is not supported");
2729 gmx_shellfc_t
*shellfc
;
2731 em_state_t state_work
{};
2733 /* Init em and store the local state in state_minimum */
2734 init_em(fplog
, NM
, cr
, outputProvider
, inputrec
,
2735 state_global
, top_global
, &state_work
, &top
,
2736 nrnb
, mu_tot
, fr
, &enerd
, &graph
, mdatoms
, &gstat
,
2737 vsite
, constr
, &shellfc
,
2738 nfile
, fnm
, &outf
, nullptr, imdport
, Flags
, wcycle
);
2740 std::vector
<size_t> atom_index
= get_atom_index(top_global
);
2741 snew(fneg
, atom_index
.size());
2742 snew(dfdx
, atom_index
.size());
2748 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2749 " which MIGHT not be accurate enough for normal mode analysis.\n"
2750 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
2751 " are fairly modest even if you recompile in double precision.\n\n");
2755 /* Check if we can/should use sparse storage format.
2757 * Sparse format is only useful when the Hessian itself is sparse, which it
2758 * will be when we use a cutoff.
2759 * For small systems (n<1000) it is easier to always use full matrix format, though.
2761 if (EEL_FULL(fr
->eeltype
) || fr
->rlist
== 0.0)
2763 GMX_LOG(mdlog
.warning
).appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2766 else if (atom_index
.size() < 1000)
2768 GMX_LOG(mdlog
.warning
).appendTextFormatted("Small system size (N=%d), using full Hessian format.",
2774 GMX_LOG(mdlog
.warning
).appendText("Using compressed symmetric sparse Hessian format.");
2778 /* Number of dimensions, based on real atoms, that is not vsites or shell */
2779 sz
= DIM
*atom_index
.size();
2781 fprintf(stderr
, "Allocating Hessian memory...\n\n");
2785 sparse_matrix
= gmx_sparsematrix_init(sz
);
2786 sparse_matrix
->compressed_symmetric
= TRUE
;
2790 snew(full_matrix
, sz
*sz
);
2797 /* Write start time and temperature */
2798 print_em_start(fplog
, cr
, walltime_accounting
, wcycle
, NM
);
2800 /* fudge nr of steps to nr of atoms */
2801 inputrec
->nsteps
= atom_index
.size()*2;
2805 fprintf(stderr
, "starting normal mode calculation '%s'\n%d steps.\n\n",
2806 *(top_global
->name
), (int)inputrec
->nsteps
);
2809 nnodes
= cr
->nnodes
;
2811 /* Make evaluate_energy do a single node force calculation */
2813 evaluate_energy(fplog
, cr
,
2814 top_global
, &state_work
, top
,
2815 inputrec
, nrnb
, wcycle
, gstat
,
2816 vsite
, constr
, fcd
, graph
, mdatoms
, fr
,
2817 mu_tot
, enerd
, vir
, pres
, -1, TRUE
);
2818 cr
->nnodes
= nnodes
;
2820 /* if forces are not small, warn user */
2821 get_state_f_norm_max(cr
, &(inputrec
->opts
), mdatoms
, &state_work
);
2823 GMX_LOG(mdlog
.warning
).appendTextFormatted("Maximum force:%12.5e", state_work
.fmax
);
2824 if (state_work
.fmax
> 1.0e-3)
2826 GMX_LOG(mdlog
.warning
).appendText(
2827 "The force is probably not small enough to "
2828 "ensure that you are at a minimum.\n"
2829 "Be aware that negative eigenvalues may occur\n"
2830 "when the resulting matrix is diagonalized.");
2833 /***********************************************************
2835 * Loop over all pairs in matrix
2837 * do_force called twice. Once with positive and
2838 * once with negative displacement
2840 ************************************************************/
2842 /* Steps are divided one by one over the nodes */
2844 for (unsigned int aid
= cr
->nodeid
; aid
< atom_index
.size(); aid
+= nnodes
)
2846 size_t atom
= atom_index
[aid
];
2847 for (size_t d
= 0; d
< DIM
; d
++)
2849 gmx_bool bBornRadii
= FALSE
;
2850 gmx_int64_t step
= 0;
2851 int force_flags
= GMX_FORCE_STATECHANGED
| GMX_FORCE_ALLFORCES
;
2854 x_min
= state_work
.s
.x
[atom
][d
];
2856 for (unsigned int dx
= 0; (dx
< 2); dx
++)
2860 state_work
.s
.x
[atom
][d
] = x_min
- der_range
;
2864 state_work
.s
.x
[atom
][d
] = x_min
+ der_range
;
2867 /* Make evaluate_energy do a single node force calculation */
2871 /* Now is the time to relax the shells */
2872 (void) relax_shell_flexcon(fplog
, cr
, bVerbose
, step
,
2873 inputrec
, bNS
, force_flags
,
2876 &state_work
.s
, &state_work
.f
, vir
, mdatoms
,
2877 nrnb
, wcycle
, graph
, &top_global
->groups
,
2878 shellfc
, fr
, bBornRadii
, t
, mu_tot
,
2885 evaluate_energy(fplog
, cr
,
2886 top_global
, &state_work
, top
,
2887 inputrec
, nrnb
, wcycle
, gstat
,
2888 vsite
, constr
, fcd
, graph
, mdatoms
, fr
,
2889 mu_tot
, enerd
, vir
, pres
, atom
*2+dx
, FALSE
);
2892 cr
->nnodes
= nnodes
;
2896 for (size_t i
= 0; i
< atom_index
.size(); i
++)
2898 copy_rvec(state_work
.f
[atom_index
[i
]], fneg
[i
]);
2903 /* x is restored to original */
2904 state_work
.s
.x
[atom
][d
] = x_min
;
2906 for (size_t j
= 0; j
< atom_index
.size(); j
++)
2908 for (size_t k
= 0; (k
< DIM
); k
++)
2911 -(state_work
.f
[atom_index
[j
]][k
] - fneg
[j
][k
])/(2*der_range
);
2918 #define mpi_type GMX_MPI_REAL
2919 MPI_Send(dfdx
[0], atom_index
.size()*DIM
, mpi_type
, MASTER(cr
),
2920 cr
->nodeid
, cr
->mpi_comm_mygroup
);
2925 for (node
= 0; (node
< nnodes
&& atom
+node
< atom_index
.size()); node
++)
2931 MPI_Recv(dfdx
[0], atom_index
.size()*DIM
, mpi_type
, node
, node
,
2932 cr
->mpi_comm_mygroup
, &stat
);
2937 row
= (atom
+ node
)*DIM
+ d
;
2939 for (size_t j
= 0; j
< atom_index
.size(); j
++)
2941 for (size_t k
= 0; k
< DIM
; k
++)
2947 if (col
>= row
&& dfdx
[j
][k
] != 0.0)
2949 gmx_sparsematrix_increment_value(sparse_matrix
,
2950 row
, col
, dfdx
[j
][k
]);
2955 full_matrix
[row
*sz
+col
] = dfdx
[j
][k
];
2962 if (bVerbose
&& fplog
)
2967 /* write progress */
2968 if (bIsMaster
&& bVerbose
)
2970 fprintf(stderr
, "\rFinished step %d out of %d",
2971 static_cast<int>(std::min(atom
+nnodes
, atom_index
.size())),
2972 static_cast<int>(atom_index
.size()));
2979 fprintf(stderr
, "\n\nWriting Hessian...\n");
2980 gmx_mtxio_write(ftp2fn(efMTX
, nfile
, fnm
), sz
, sz
, full_matrix
, sparse_matrix
);
2983 finish_em(cr
, outf
, walltime_accounting
, wcycle
);
2985 walltime_accounting_set_nsteps_done(walltime_accounting
, atom_index
.size()*2);