2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016, 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.
37 * Implements routines in boxutilities.h.
39 * Utility functions for handling boxes.
43 #include "boxutilities.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/math/vectypes.h"
51 #include "gromacs/mdtypes/inputrec.h"
52 #include "gromacs/mdtypes/md_enums.h"
53 #include "gromacs/mdtypes/state.h"
54 #include "gromacs/pbcutil/pbc.h"
55 #include "gromacs/utility/basedefinitions.h"
56 #include "gromacs/utility/real.h"
58 /*! \brief Change box components to preserve the relative box shape
60 * Change box components to box[XX][XX]*box_rel to preserve the relative box shape
62 static void do_box_rel(const t_inputrec *ir, matrix box_rel,
63 matrix b, gmx_bool bInit)
67 for (d = YY; d <= ZZ; d++)
69 for (d2 = XX; d2 <= (ir->epct == epctSEMIISOTROPIC ? YY : ZZ); d2++)
71 /* We need to check if this box component is deformed
72 * or if deformation of another component might cause
73 * changes in this component due to box corrections.
75 if (ir->deform[d][d2] == 0 &&
76 !(d == ZZ && d2 == XX && ir->deform[d][YY] != 0 &&
77 (b[YY][d2] != 0 || ir->deform[YY][d2] != 0)))
81 box_rel[d][d2] = b[d][d2]/b[XX][XX];
85 b[d][d2] = b[XX][XX]*box_rel[d][d2];
92 void preserve_box_shape(const t_inputrec *ir, matrix box_rel, matrix b)
94 if (inputrecPreserveShape(ir))
96 do_box_rel(ir, box_rel, b, FALSE);
100 void set_box_rel(const t_inputrec *ir, t_state *state)
102 /* Make sure the box obeys the restrictions before we fix the ratios */
103 correct_box(NULL, 0, state->box, NULL);
105 clear_mat(state->box_rel);
107 if (inputrecPreserveShape(ir))
109 do_box_rel(ir, state->box_rel, state->box, TRUE);