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.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 /* This file is completely threadsafe - keep it that way! */
47 #include "gromacs/math/paddedvector.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/math/veccompare.h"
50 #include "gromacs/mdtypes/awh_history.h"
51 #include "gromacs/mdtypes/df_history.h"
52 #include "gromacs/mdtypes/inputrec.h"
53 #include "gromacs/mdtypes/md_enums.h"
54 #include "gromacs/mdtypes/pull_params.h"
55 #include "gromacs/mdtypes/swaphistory.h"
56 #include "gromacs/pbcutil/boxutilities.h"
57 #include "gromacs/pbcutil/pbc.h"
58 #include "gromacs/utility/compare.h"
59 #include "gromacs/utility/gmxassert.h"
60 #include "gromacs/utility/smalloc.h"
62 /* The source code in this file should be thread-safe.
63 Please keep it that way. */
65 history_t::history_t() :
68 disre_rm3tav(nullptr),
71 orire_Dtav(nullptr){};
73 ekinstate_t::ekinstate_t() :
84 clear_mat(ekin_total);
87 void init_gtc_state(t_state* state, int ngtc, int nnhpres, int nhchainlength)
90 state->nnhpres = nnhpres;
91 state->nhchainlength = nhchainlength;
92 state->nosehoover_xi.resize(state->nhchainlength * state->ngtc, 0);
93 state->nosehoover_vxi.resize(state->nhchainlength * state->ngtc, 0);
94 state->therm_integral.resize(state->ngtc, 0);
95 state->baros_integral = 0.0;
96 state->nhpres_xi.resize(state->nhchainlength * nnhpres, 0);
97 state->nhpres_vxi.resize(state->nhchainlength * nnhpres, 0);
101 /* Checkpoint code relies on this function having no effect if
102 state->natoms is > 0 and passed as natoms. */
103 void state_change_natoms(t_state* state, int natoms)
105 state->natoms = natoms;
107 /* We need padding, since we might use SIMD access, but the
108 * containers here all ensure that. */
109 if (state->flags & (1 << estX))
111 state->x.resizeWithPadding(natoms);
113 if (state->flags & (1 << estV))
115 state->v.resizeWithPadding(natoms);
117 if (state->flags & (1 << estCGP))
119 state->cg_p.resizeWithPadding(natoms);
123 void init_dfhist_state(t_state* state, int dfhistNumLambda)
125 if (dfhistNumLambda > 0)
127 snew(state->dfhist, 1);
128 init_df_history(state->dfhist, dfhistNumLambda);
132 state->dfhist = nullptr;
136 void comp_state(const t_state* st1, const t_state* st2, gmx_bool bRMSD, real ftol, real abstol)
140 fprintf(stdout, "comparing flags\n");
141 cmp_int(stdout, "flags", -1, st1->flags, st2->flags);
142 fprintf(stdout, "comparing box\n");
143 cmp_rvecs(stdout, "box", DIM, st1->box, st2->box, FALSE, ftol, abstol);
144 fprintf(stdout, "comparing box_rel\n");
145 cmp_rvecs(stdout, "box_rel", DIM, st1->box_rel, st2->box_rel, FALSE, ftol, abstol);
146 fprintf(stdout, "comparing boxv\n");
147 cmp_rvecs(stdout, "boxv", DIM, st1->boxv, st2->boxv, FALSE, ftol, abstol);
148 if (st1->flags & (1 << estSVIR_PREV))
150 fprintf(stdout, "comparing shake vir_prev\n");
151 cmp_rvecs(stdout, "svir_prev", DIM, st1->svir_prev, st2->svir_prev, FALSE, ftol, abstol);
153 if (st1->flags & (1 << estFVIR_PREV))
155 fprintf(stdout, "comparing force vir_prev\n");
156 cmp_rvecs(stdout, "fvir_prev", DIM, st1->fvir_prev, st2->fvir_prev, FALSE, ftol, abstol);
158 if (st1->flags & (1 << estPRES_PREV))
160 fprintf(stdout, "comparing prev_pres\n");
161 cmp_rvecs(stdout, "pres_prev", DIM, st1->pres_prev, st2->pres_prev, FALSE, ftol, abstol);
163 cmp_int(stdout, "ngtc", -1, st1->ngtc, st2->ngtc);
164 cmp_int(stdout, "nhchainlength", -1, st1->nhchainlength, st2->nhchainlength);
165 if (st1->ngtc == st2->ngtc && st1->nhchainlength == st2->nhchainlength)
167 for (i = 0; i < st1->ngtc; i++)
169 nc = i * st1->nhchainlength;
170 for (j = 0; j < nc; j++)
172 cmp_real(stdout, "nosehoover_xi", i, st1->nosehoover_xi[nc + j],
173 st2->nosehoover_xi[nc + j], ftol, abstol);
177 cmp_int(stdout, "nnhpres", -1, st1->nnhpres, st2->nnhpres);
178 if (st1->nnhpres == st2->nnhpres && st1->nhchainlength == st2->nhchainlength)
180 for (i = 0; i < st1->nnhpres; i++)
182 nc = i * st1->nhchainlength;
183 for (j = 0; j < nc; j++)
185 cmp_real(stdout, "nosehoover_xi", i, st1->nhpres_xi[nc + j], st2->nhpres_xi[nc + j],
191 cmp_int(stdout, "natoms", -1, st1->natoms, st2->natoms);
192 if (st1->natoms == st2->natoms)
194 if ((st1->flags & (1 << estX)) && (st2->flags & (1 << estX)))
196 fprintf(stdout, "comparing x\n");
197 cmp_rvecs(stdout, "x", st1->natoms, st1->x.rvec_array(), st2->x.rvec_array(), bRMSD,
200 if ((st1->flags & (1 << estV)) && (st2->flags & (1 << estV)))
202 fprintf(stdout, "comparing v\n");
203 cmp_rvecs(stdout, "v", st1->natoms, st1->v.rvec_array(), st2->v.rvec_array(), bRMSD,
209 rvec* makeRvecArray(gmx::ArrayRef<const gmx::RVec> v, gmx::index n)
211 GMX_ASSERT(v.ssize() >= n, "We can't copy more elements than the vector size");
217 const rvec* vPtr = as_rvec_array(v.data());
218 for (gmx::index i = 0; i < n; i++)
220 copy_rvec(vPtr[i], dest[i]);
247 // It would be nicer to initialize these with {} or {{0}} in the
248 // above initialization list, but uncrustify doesn't understand
250 // TODO Fix this if we switch to clang-format some time.
255 clear_mat(pres_prev);
256 clear_mat(svir_prev);
257 clear_mat(fvir_prev);
260 void set_box_rel(const t_inputrec* ir, t_state* state)
262 /* Make sure the box obeys the restrictions before we fix the ratios */
263 correct_box(nullptr, 0, state->box);
265 clear_mat(state->box_rel);
267 if (inputrecPreserveShape(ir))
269 const int ndim = ir->epct == epctSEMIISOTROPIC ? 2 : 3;
270 do_box_rel(ndim, ir->deform, state->box_rel, state->box, true);
274 void preserve_box_shape(const t_inputrec* ir, matrix box_rel, matrix box)
276 if (inputrecPreserveShape(ir))
278 const int ndim = ir->epct == epctSEMIISOTROPIC ? 2 : 3;
279 do_box_rel(ndim, ir->deform, box_rel, box, false);
283 void initialize_lambdas(FILE* fplog,
284 const t_inputrec& ir,
287 gmx::ArrayRef<real> lambda,
290 /* TODO: Clean up initialization of fep_state and lambda in
291 t_state. This function works, but could probably use a logic
292 rewrite to keep all the different types of efep straight. */
294 if ((ir.efep == efepNO) && (!ir.bSimTemp))
299 const t_lambda* fep = ir.fepvals;
302 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
303 if checkpoint is set -- a kludge is in for now
307 for (int i = 0; i < efptNR; i++)
310 /* overwrite lambda state with init_lambda for now for backwards compatibility */
311 if (fep->init_lambda >= 0) /* if it's -1, it was never initialized */
313 thisLambda = fep->init_lambda;
317 thisLambda = fep->all_lambda[i][fep->init_fep_state];
321 lambda[i] = thisLambda;
325 lam0[i] = thisLambda;
330 /* need to rescale control temperatures to match current state */
331 for (int i = 0; i < ir.opts.ngtc; i++)
333 if (ir.opts.ref_t[i] > 0)
335 ir.opts.ref_t[i] = ir.simtempvals->temperatures[fep->init_fep_state];
340 /* Send to the log the information on the current lambdas */
341 if (fplog != nullptr)
343 fprintf(fplog, "Initial vector of lambda components:[ ");
344 for (const auto& l : lambda)
346 fprintf(fplog, "%10.4f ", l);
348 fprintf(fplog, "]\n");