Added basic tests for gmx msd
[alexxy/gromacs.git] / src / gromacs / gmxana / tests / gmx_msd.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
37  * Tests for gmx msd.
38  *
39  * \author Kevin Boyd <kevin.boyd@uconn.edu>
40  */
41
42 #include "gmxpre.h"
43
44 #include "gromacs/gmxana/gmx_ana.h"
45 #include "gromacs/utility/futil.h"
46 #include "gromacs/utility/textreader.h"
47
48 #include "testutils/cmdlinetest.h"
49 #include "testutils/refdata.h"
50 #include "testutils/testfilemanager.h"
51 #include "testutils/textblockmatchers.h"
52 #include "testutils/xvgtest.h"
53
54 namespace
55 {
56
57 using gmx::test::CommandLine;
58 using gmx::test::XvgMatch;
59
60 class MsdTest : public gmx::test::CommandLineTestBase
61 {
62     public:
63         MsdTest()
64         {
65             setOutputFile("-o", "msd.xvg", XvgMatch());
66             setInputFile("-f", "msd_traj.xtc");
67             setInputFile("-s", "msd_coords.gro");
68             setInputFile("-n", "msd.ndx");
69         }
70
71         void runTest(const CommandLine &args)
72         {
73             CommandLine &cmdline = commandLine();
74             cmdline.merge(args);
75             ASSERT_EQ(0, gmx_msd(cmdline.argc(), cmdline.argv()));
76             checkOutputFiles();
77         }
78 };
79
80 /* msd_traj.xtc contains a 10 frame (1 ps per frame) simulation
81  * containing 3 atoms, with different starting positions but identical
82  * displacements. The displacements are calculated to yield the following
83  * diffusion coefficients when lag is calculated ONLY FROM TIME 0
84  * D_x = 8 * 10 ^ -5 cm^2 /s, D_y = 4 * 10^ -5 cm^2 /s , D_z = 0
85  *
86  * To test for these results, -trestart is set to a larger value than the
87  * total simulation length, so that only lag 0 is calculated
88  */
89
90 // for 3D, (8 + 4 + 0) / 3 should yield 4 cm^2 / s
91 TEST_F(MsdTest, threeDimensionalDiffusion)
92 {
93     const char *const cmdline[] = {
94         "msd", "-mw", "no", "-trestart", "200",
95     };
96     runTest(CommandLine(cmdline));
97 }
98
99 // for lateral z, (8 + 4) / 2 should yield 6 cm^2 /s
100 TEST_F(MsdTest, twoDimensionalDiffusion)
101 {
102     const char *const cmdline[] = {
103         "msd", "-mw", "no", "-trestart", "200", "-lateral", "z"
104     };
105     runTest(CommandLine(cmdline));
106 }
107
108 // for type x, should yield 8 cm^2 / s
109 TEST_F(MsdTest, oneDimensionalDiffusion)
110 {
111     const char *const cmdline[] = {
112         "msd", "-mw", "no", "-trestart", "200", "-type", "x"
113     };
114     runTest(CommandLine(cmdline));
115 }
116 } //namespace