Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / pbcutil / boxutilities.cpp
blob4ee39267ce6f00c2cea591cb84c3c2fa34f886a3
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2017,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.
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.
35 /*! \internal \file
36 * \brief
37 * Implements routines in boxutilities.h.
39 * Utility functions for handling boxes.
41 #include "gmxpre.h"
43 #include "boxutilities.h"
45 #include <algorithm>
47 #include "gromacs/math/utilities.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/math/vectypes.h"
51 void do_box_rel(int ndim, const matrix deform, matrix box_rel,
52 matrix b, bool bInit)
54 for (int d = YY; d <= ZZ; ++d)
56 for (int d2 = XX; d2 < ndim; ++d2)
58 /* We need to check if this box component is deformed
59 * or if deformation of another component might cause
60 * changes in this component due to box corrections.
62 if (deform[d][d2] == 0 &&
63 !(d == ZZ && d2 == XX && deform[d][YY] != 0 &&
64 (b[YY][d2] != 0 || deform[YY][d2] != 0)))
66 if (bInit)
68 box_rel[d][d2] = b[d][d2]/b[XX][XX];
70 else
72 b[d][d2] = b[XX][XX]*box_rel[d][d2];
79 namespace gmx
82 namespace
85 //! Whether two box elements are equal (with a tolerance).
86 bool boxElementEqual(real element1, real element2)
88 // Compare with a relative tolerance (for big boxes) and with
89 // an absolute tolerance (small boxes are generally not specified with very
90 // high number of decimals).
91 return gmx_within_tol(element1, element2, 10*GMX_REAL_EPS)
92 || std::fabs(element1 - element2) < 1e-3;
95 } // namespace
97 bool boxesAreEqual(const matrix box1, const matrix box2)
99 return boxElementEqual(box1[XX][XX], box2[XX][XX])
100 && boxElementEqual(box1[YY][XX], box2[YY][XX])
101 && boxElementEqual(box1[YY][YY], box2[YY][YY])
102 && boxElementEqual(box1[ZZ][XX], box2[ZZ][XX])
103 && boxElementEqual(box1[ZZ][YY], box2[ZZ][YY])
104 && boxElementEqual(box1[ZZ][ZZ], box2[ZZ][ZZ]);
107 bool boxIsZero(const matrix box)
109 return boxElementEqual(box[XX][XX], 0.0)
110 && boxElementEqual(box[YY][XX], 0.0)
111 && boxElementEqual(box[YY][YY], 0.0)
112 && boxElementEqual(box[ZZ][XX], 0.0)
113 && boxElementEqual(box[ZZ][YY], 0.0)
114 && boxElementEqual(box[ZZ][ZZ], 0.0);
117 } // namespace gmx