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,2018, 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.
37 /* This file is completely threadsafe - keep it that way! */
46 #include "gromacs/math/paddedvector.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/math/veccompare.h"
49 #include "gromacs/mdtypes/awh-history.h"
50 #include "gromacs/mdtypes/df_history.h"
51 #include "gromacs/mdtypes/inputrec.h"
52 #include "gromacs/mdtypes/md_enums.h"
53 #include "gromacs/mdtypes/pull-params.h"
54 #include "gromacs/mdtypes/swaphistory.h"
55 #include "gromacs/pbcutil/boxutilities.h"
56 #include "gromacs/pbcutil/pbc.h"
57 #include "gromacs/utility/compare.h"
58 #include "gromacs/utility/gmxassert.h"
59 #include "gromacs/utility/smalloc.h"
61 /* The source code in this file should be thread-safe.
62 Please keep it that way. */
64 history_t::history_t() : disre_initf(0),
66 disre_rm3tav(nullptr),
73 ekinstate_t::ekinstate_t() : bUpToDate(FALSE),
83 clear_mat(ekin_total);
86 void init_gtc_state(t_state *state, int ngtc, int nnhpres, int nhchainlength)
89 state->nnhpres = nnhpres;
90 state->nhchainlength = nhchainlength;
91 state->nosehoover_xi.resize(state->nhchainlength*state->ngtc, 0);
92 state->nosehoover_vxi.resize(state->nhchainlength*state->ngtc, 0);
93 state->therm_integral.resize(state->ngtc, 0);
94 state->baros_integral = 0.0;
95 state->nhpres_xi.resize(state->nhchainlength*nnhpres, 0);
96 state->nhpres_vxi.resize(state->nhchainlength*nnhpres, 0);
100 /* Checkpoint code relies on this function having no effect if
101 state->natoms is > 0 and passed as natoms. */
102 void state_change_natoms(t_state *state, int natoms)
104 state->natoms = natoms;
106 /* We need padding, since we might use SIMD access, but the
107 * containers here all ensure that. */
108 if (state->flags & (1 << estX))
110 state->x.resizeWithPadding(natoms);
112 if (state->flags & (1 << estV))
114 state->v.resizeWithPadding(natoms);
116 if (state->flags & (1 << estCGP))
118 state->cg_p.resizeWithPadding(natoms);
122 void init_dfhist_state(t_state *state, int dfhistNumLambda)
124 if (dfhistNumLambda > 0)
126 snew(state->dfhist, 1);
127 init_df_history(state->dfhist, dfhistNumLambda);
131 state->dfhist = nullptr;
135 void comp_state(const t_state *st1, const t_state *st2,
136 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",
173 i, st1->nosehoover_xi[nc+j], 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",
186 i, st1->nhpres_xi[nc+j], st2->nhpres_xi[nc+j], ftol, abstol);
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, ftol, abstol);
199 if ((st1->flags & (1<<estV)) && (st2->flags & (1<<estV)))
201 fprintf(stdout, "comparing v\n");
202 cmp_rvecs(stdout, "v", st1->natoms, st1->v.rvec_array(), st2->v.rvec_array(), bRMSD, ftol, abstol);
207 rvec *makeRvecArray(gmx::ArrayRef<const gmx::RVec> v,
210 GMX_ASSERT(v.size() >= n, "We can't copy more elements than the vector size");
216 const rvec *vPtr = as_rvec_array(v.data());
217 for (gmx::index i = 0; i < n; i++)
219 copy_rvec(vPtr[i], dest[i]);
225 t_state::t_state() : natoms(0),
245 // It would be nicer to initialize these with {} or {{0}} in the
246 // above initialization list, but uncrustify doesn't understand
248 // TODO Fix this if we switch to clang-format some time.
253 clear_mat(pres_prev);
254 clear_mat(svir_prev);
255 clear_mat(fvir_prev);
258 void set_box_rel(const t_inputrec *ir, t_state *state)
260 /* Make sure the box obeys the restrictions before we fix the ratios */
261 correct_box(nullptr, 0, state->box, nullptr);
263 clear_mat(state->box_rel);
265 if (inputrecPreserveShape(ir))
267 const int ndim = ir->epct == epctSEMIISOTROPIC ? 2 : 3;
268 do_box_rel(ndim, ir->deform, state->box_rel, state->box, true);
272 void preserve_box_shape(const t_inputrec *ir, matrix box_rel, matrix box)
274 if (inputrecPreserveShape(ir))
276 const int ndim = ir->epct == epctSEMIISOTROPIC ? 2 : 3;
277 do_box_rel(ndim, ir->deform, box_rel, box, false);