2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 * \brief Defines the box deformation code.
38 * \todo The .mdp specification should have a boolean for this module.
40 * \todo grompp should set up fields in the tpr file that carry the
41 * information about the original box, then the deform module would
42 * can be build alongside the update module, rather than need to be
43 * set up before the checkpoint is read.
45 * \author Berk Hess <hess@kth.se>
46 * \author Mark Abraham <mark.j.abraham@gmail.com>
47 * \ingroup module_mdlib
51 #include "boxdeformation.h"
55 #include "gromacs/gmxlib/network.h"
56 #include "gromacs/math/invertmatrix.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/mdtypes/commrec.h"
59 #include "gromacs/mdtypes/inputrec.h"
60 #include "gromacs/mdtypes/md_enums.h"
61 #include "gromacs/utility/exceptions.h"
66 std::unique_ptr<BoxDeformation> prepareBoxDeformation(const matrix& initialBox,
68 const t_inputrec& inputrec)
70 if (!inputrecDeform(&inputrec))
74 if (!EI_DYNAMICS(inputrec.eI))
76 GMX_THROW(NotImplementedError(
77 "Box deformation is only supported with dynamical integrators"));
81 // Only the rank that read the tpr has the global state, and thus
82 // the initial box, so we pass that around.
85 copy_mat(initialBox, box);
89 gmx_bcast(sizeof(box), box, cr);
92 return std::make_unique<BoxDeformation>(inputrec.delta_t, inputrec.init_step, inputrec.deform, box);
95 BoxDeformation::BoxDeformation(double timeStep,
97 const tensor& deformationTensor,
98 const matrix& referenceBox) :
100 initialStep_(initialStep)
102 copy_mat(deformationTensor, deformationTensor_);
103 copy_mat(referenceBox, referenceBox_);
106 void BoxDeformation::apply(ArrayRef<RVec> x, matrix box, int64_t step)
108 matrix updatedBox, invbox, mu;
110 double elapsedTime = (step + 1 - initialStep_) * timeStep_;
111 copy_mat(box, updatedBox);
112 for (int i = 0; i < DIM; i++)
114 for (int j = 0; j < DIM; j++)
116 if (deformationTensor_[i][j] != 0)
118 updatedBox[i][j] = referenceBox_[i][j] + elapsedTime * deformationTensor_[i][j];
122 /* We correct the off-diagonal elements,
123 * which can grow indefinitely during shearing,
124 * so the shifts do not get messed up.
126 for (int i = 1; i < DIM; i++)
128 for (int j = i - 1; j >= 0; j--)
130 while (updatedBox[i][j] - box[i][j] > 0.5 * updatedBox[j][j])
132 rvec_dec(updatedBox[i], updatedBox[j]);
134 while (updatedBox[i][j] - box[i][j] < -0.5 * updatedBox[j][j])
136 rvec_inc(updatedBox[i], updatedBox[j]);
140 invertBoxMatrix(box, invbox);
141 // Return the updated box
142 copy_mat(updatedBox, box);
143 mmul_ur0(box, invbox, mu);
145 for (auto& thisX : x)
147 thisX[XX] = mu[XX][XX] * thisX[XX] + mu[YY][XX] * thisX[YY] + mu[ZZ][XX] * thisX[ZZ];
148 thisX[YY] = mu[YY][YY] * thisX[YY] + mu[ZZ][YY] * thisX[ZZ];
149 thisX[ZZ] = mu[ZZ][ZZ] * thisX[ZZ];