Remove gmx custom fixed int (e.g. gmx_int64_t) types
[alexxy/gromacs.git] / src / gromacs / mdlib / boxdeformation.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  * \brief Defines the box deformation code.
37  *
38  * \todo The .mdp specification should have a boolean for this module.
39  *
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.
44  *
45  * \author Berk Hess <hess@kth.se>
46  * \author Mark Abraham <mark.j.abraham@gmail.com>
47  * \ingroup module_mdlib
48  */
49 #include "gmxpre.h"
50
51 #include "boxdeformation.h"
52
53 #include "gromacs/compat/make_unique.h"
54 #include "gromacs/gmxlib/network.h"
55 #include "gromacs/math/invertmatrix.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/mdtypes/commrec.h"
58 #include "gromacs/mdtypes/inputrec.h"
59 #include "gromacs/mdtypes/md_enums.h"
60 #include "gromacs/utility/exceptions.h"
61
62 namespace gmx
63 {
64
65 std::unique_ptr<BoxDeformation>
66 prepareBoxDeformation(const matrix     &initialBox,
67                       t_commrec        *cr,
68                       const t_inputrec &inputrec)
69 {
70     if (!inputrecDeform(&inputrec))
71     {
72         return nullptr;
73     }
74     if (!EI_DYNAMICS(inputrec.eI))
75     {
76         GMX_THROW(NotImplementedError("Box deformation is only supported with dynamical integrators"));
77     }
78
79     matrix box;
80     // Only the rank that read the tpr has the global state, and thus
81     // the initial box, so we pass that around.
82     if (SIMMASTER(cr))
83     {
84         copy_mat(initialBox, box);
85     }
86     if (PAR(cr))
87     {
88         gmx_bcast(sizeof(box), box, cr);
89     }
90
91     return compat::make_unique<BoxDeformation>(inputrec.delta_t,
92                                                inputrec.init_step,
93                                                inputrec.deform,
94                                                box);
95 }
96
97 BoxDeformation::BoxDeformation(double        timeStep,
98                                int64_t       initialStep,
99                                const tensor &deformationTensor,
100                                const matrix &referenceBox)
101     : timeStep_(timeStep),
102       initialStep_(initialStep)
103 {
104     copy_mat(deformationTensor, deformationTensor_);
105     copy_mat(referenceBox, referenceBox_);
106 }
107
108 void
109 BoxDeformation::apply(ArrayRef<RVec> x,
110                       matrix         box,
111                       int64_t        step)
112 {
113     matrix updatedBox, invbox, mu;
114
115     double elapsedTime = (step + 1 - initialStep_) * timeStep_;
116     copy_mat(box, updatedBox);
117     for (int i = 0; i < DIM; i++)
118     {
119         for (int j = 0; j < DIM; j++)
120         {
121             if (deformationTensor_[i][j] != 0)
122             {
123                 updatedBox[i][j] =
124                     referenceBox_[i][j] + elapsedTime * deformationTensor_[i][j];
125             }
126         }
127     }
128     /* We correct the off-diagonal elements,
129      * which can grow indefinitely during shearing,
130      * so the shifts do not get messed up.
131      */
132     for (int i = 1; i < DIM; i++)
133     {
134         for (int j = i-1; j >= 0; j--)
135         {
136             while (updatedBox[i][j] - box[i][j] > 0.5 * updatedBox[j][j])
137             {
138                 rvec_dec(updatedBox[i], updatedBox[j]);
139             }
140             while (updatedBox[i][j] - box[i][j] < -0.5 * updatedBox[j][j])
141             {
142                 rvec_inc(updatedBox[i], updatedBox[j]);
143             }
144         }
145     }
146     invertBoxMatrix(box, invbox);
147     // Return the updated box
148     copy_mat(updatedBox, box);
149     mmul_ur0(box, invbox, mu);
150
151     for (auto &thisX : x)
152     {
153         thisX[XX] = mu[XX][XX]*thisX[XX] + mu[YY][XX]*thisX[YY] + mu[ZZ][XX]*thisX[ZZ];
154         thisX[YY] = mu[YY][YY]*thisX[YY] + mu[ZZ][YY]*thisX[ZZ];
155         thisX[ZZ] = mu[ZZ][ZZ]*thisX[ZZ];
156     }
157 }
158
159 } // namespace