f32c0c0fb1a5fb5752c8bd9abde390e0c527e855
[alexxy/gromacs.git] / src / gromacs / fileio / tests / readinp.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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.
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 utilities for routines that parse fields e.g. from grompp input
38  *
39  * \author Mark Abraham <mark.j.abraham@gmail.com>
40  */
41 #include "gmxpre.h"
42
43 #include "gromacs/fileio/readinp.h"
44
45 #include <gtest/gtest.h>
46
47 #include "gromacs/fileio/warninp.h"
48 #include "gromacs/utility/smalloc.h"
49 #include "gromacs/utility/unique_cptr.h"
50
51 namespace gmx
52 {
53 namespace testing
54 {
55
56 class ReadTest : public ::testing::Test
57 {
58     public:
59         ReadTest() : numInputs_(1),
60                      inputField_(0),
61                      inpGuard_(),
62                      wi_(),
63                      wiGuard_()
64         {
65             snew(inputField_, numInputs_);
66             inpGuard_.reset(inputField_);
67
68             inputField_[0].count     = 0;
69             inputField_[0].bObsolete = FALSE;
70             inputField_[0].bSet      = FALSE;
71             inputField_[0].name      = (char *) "test";
72             inputField_[0].inp_count = 0;
73
74             wi_ = init_warning(FALSE, 0);
75             wiGuard_.reset(wi_);
76         }
77
78         int                                            numInputs_;
79         t_inpfile                                     *inputField_;
80         gmx::unique_cptr<t_inpfile>                    inpGuard_;
81         warninp_t                                      wi_;
82         gmx::unique_cptr<struct warninp, free_warning> wiGuard_;
83 };
84
85 TEST_F(ReadTest, get_eint_ReadsInteger)
86 {
87     inputField_[0].value = (char *) "1";
88     ASSERT_EQ(1, get_eint(&numInputs_, &inputField_, "test", 2, wi_));
89     ASSERT_FALSE(warning_errors_exist(wi_));
90 }
91
92 TEST_F(ReadTest, get_eint_WarnsAboutFloat)
93 {
94     inputField_[0].value = (char *) "0.8";
95     get_eint(&numInputs_, &inputField_, "test", 2, wi_);
96     ASSERT_TRUE(warning_errors_exist(wi_));
97 }
98
99 TEST_F(ReadTest, get_eint_WarnsAboutString)
100 {
101     inputField_[0].value = (char *) "hello";
102     get_eint(&numInputs_, &inputField_, "test", 2, wi_);
103     ASSERT_TRUE(warning_errors_exist(wi_));
104 }
105
106 TEST_F(ReadTest, get_eint64_ReadsInteger)
107 {
108     inputField_[0].value = (char *) "1";
109     ASSERT_EQ(1, get_eint64(&numInputs_, &inputField_, "test", 2, wi_));
110     ASSERT_FALSE(warning_errors_exist(wi_));
111 }
112
113 TEST_F(ReadTest, get_eint64_WarnsAboutFloat)
114 {
115     inputField_[0].value = (char *) "0.8";
116     get_eint64(&numInputs_, &inputField_, "test", 2, wi_);
117     ASSERT_TRUE(warning_errors_exist(wi_));
118 }
119
120 TEST_F(ReadTest, get_eint64_WarnsAboutString)
121 {
122     inputField_[0].value = (char *) "hello";
123     get_eint64(&numInputs_, &inputField_, "test", 2, wi_);
124     ASSERT_TRUE(warning_errors_exist(wi_));
125 }
126
127 TEST_F(ReadTest, get_ereal_ReadsInteger)
128 {
129     inputField_[0].value = (char *) "1";
130     ASSERT_EQ(1, get_ereal(&numInputs_, &inputField_, "test", 2, wi_));
131     ASSERT_FALSE(warning_errors_exist(wi_));
132 }
133
134 TEST_F(ReadTest, get_ereal_ReadsFloat)
135 {
136     inputField_[0].value = (char *) "0.8";
137     ASSERT_EQ(0.8, get_ereal(&numInputs_, &inputField_, "test", 2, wi_));
138     ASSERT_FALSE(warning_errors_exist(wi_));
139 }
140
141 TEST_F(ReadTest, get_ereal_WarnsAboutString)
142 {
143     inputField_[0].value = (char *) "hello";
144     get_ereal(&numInputs_, &inputField_, "test", 2, wi_);
145     ASSERT_TRUE(warning_errors_exist(wi_));
146 }
147
148 } // namespace
149 } // namespace