5d5bd967df52f86794f1275f8360863d8d616f7c
[alexxy/gromacs.git] / src / gromacs / correlationfunctions / tests / expfit.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014,2015,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.
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  * Implements test of exponential fitting routines
38  *
39  * \author Anders G&auml;rden&auml;s <anders.gardenas@gmail.com>
40  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
41  * \ingroup module_correlationfunctions
42  */
43 #include "gmxpre.h"
44
45 #include "gromacs/correlationfunctions/expfit.h"
46
47 #include "config.h"
48
49 #include <cmath>
50
51 #include <gtest/gtest.h>
52
53 #include "gromacs/fileio/oenv.h"
54 #include "gromacs/fileio/xvgr.h"
55 #include "gromacs/utility/smalloc.h"
56
57 #include "testutils/refdata.h"
58 #include "testutils/testasserts.h"
59 #include "testutils/testfilemanager.h"
60
61 namespace gmx
62 {
63
64 namespace
65 {
66
67 class ExpfitData
68 {
69     public:
70         int               nrLines_;
71         std::vector<real> x_, y_;
72         real              startTime_, endTime_, dt_;
73 };
74
75 class ExpfitTest : public ::testing::Test
76 {
77
78     protected:
79         static std::vector<ExpfitData> data_;
80         test::TestReferenceData        refData_;
81         test::TestReferenceChecker     checker_;
82         ExpfitTest( )
83             : checker_(refData_.rootChecker())
84         {
85         }
86
87         // Static initiation, only run once every test.
88         static void SetUpTestCase()
89         {
90             double                ** tempValues = nullptr;
91             std::vector<std::string> fileName;
92             fileName.push_back(test::TestFileManager::getInputFilePath("testINVEXP.xvg"));
93             fileName.push_back(test::TestFileManager::getInputFilePath("testPRES.xvg"));
94             fileName.push_back(test::TestFileManager::getInputFilePath("testINVEXP79.xvg"));
95             fileName.push_back(test::TestFileManager::getInputFilePath("testERF.xvg"));
96             fileName.push_back(test::TestFileManager::getInputFilePath("testERREST.xvg"));
97             for (std::vector<std::string>::iterator i = fileName.begin(); i < fileName.end(); ++i)
98             {
99                 const char * name = i->c_str();
100                 int          nrColumns;
101                 ExpfitData   ed;
102                 ed.nrLines_   = read_xvg(name, &tempValues, &nrColumns);
103                 ed.dt_        = tempValues[0][1] - tempValues[0][0];
104                 ed.startTime_ = tempValues[0][0];
105                 ed.endTime_   = tempValues[0][ed.nrLines_-1];
106                 for (int j = 0; j  < ed.nrLines_; j++)
107                 {
108                     ed.x_.push_back((real)tempValues[0][j]);
109                     ed.y_.push_back((real)tempValues[1][j]);
110                 }
111                 data_.push_back(ed);
112
113                 // Free memory that was allocated in read_xvg
114                 for (int j = 0; j < nrColumns; j++)
115                 {
116                     sfree(tempValues[j]);
117                     tempValues[j] = nullptr;
118                 }
119                 sfree(tempValues);
120                 tempValues = nullptr;
121             }
122         }
123
124         static void TearDownTestCase()
125         {
126         }
127
128         void test(int type, double result[], double tolerance,
129                   unsigned int testType)
130         {
131             int               nfitparm = effnNparams(type);
132             gmx_output_env_t *oenv;
133
134             if (testType >= data_.size())
135             {
136                 GMX_THROW(InvalidInputError("testType out of range"));
137             }
138             output_env_init_default(&oenv);
139             do_lmfit(data_[testType].nrLines_,
140                      &(data_[testType].y_[0]),
141                      nullptr,
142                      data_[testType].dt_,
143                      &(data_[testType].x_[0]),
144                      data_[testType].startTime_,
145                      data_[testType].endTime_,
146                      oenv, false, type, result, 0, nullptr);
147             output_env_done(oenv);
148             checker_.setDefaultTolerance(test::relativeToleranceAsFloatingPoint(1, tolerance));
149             checker_.checkSequenceArray(nfitparm, result, "result");
150         }
151 };
152
153
154 //static var
155 std::vector<ExpfitData> ExpfitTest::data_;
156
157 // TODO calling test() leads to a fatal error, which we could in
158 // principle test for.
159 #if HAVE_LMFIT
160
161 TEST_F (ExpfitTest, EffnEXP1) {
162     double  param[] = {25};
163     test(effnEXP1, param, 1e-5, 0);
164 }
165
166 TEST_F (ExpfitTest, EffnEXP2) {
167     double  param[] = {35, 0.5};
168     test(effnEXP2, param, 3e-5, 0);
169 }
170
171 TEST_F (ExpfitTest, EffnEXPEXP) {
172     double param[] = {5, 0.5, 45};
173     test(effnEXPEXP, param, 1e-2, 0);
174 }
175
176 TEST_F (ExpfitTest, EffnEXP5) {
177     double  param[] = {0.5, 5, 0.5, 50, 0.002};
178     test(effnEXP5, param, 1e-2, 2);
179 }
180
181 TEST_F (ExpfitTest, EffnEXP7) {
182     double  param[] = {0.1, 2, 0.5, 30, 0.3, 50, -0.002};
183     test(effnEXP7, param, 1e-2, 2);
184 }
185
186 TEST_F (ExpfitTest, EffnEXP9) {
187     double  param[] = {0.4, 5, 0.2, 30, 0.1, 70, 0.2, 200, -0.05};
188     test(effnEXP9, param, 4e-2, 2);
189 }
190
191 TEST_F (ExpfitTest, EffnERF) {
192     double  param[] = {80, 120, 180, 5};
193     test(effnERF, param, 1e-1, 3);
194 }
195
196 TEST_F (ExpfitTest, EffnERREST) {
197     double  param[] = {1, 0.9, 100};
198     test(effnERREST, param, 5e-3, 4);
199 }
200
201 TEST_F (ExpfitTest, EffnVAC) {
202     double param[] = {0.6, 0.1};
203     test(effnVAC, param, 0.05, 0);
204 }
205
206 TEST_F (ExpfitTest, EffnPRES) {
207     double param[] = {0.6, 10, 7, 1, 0.25, 2};
208     test(effnPRES, param, 1e-4, 1);
209 }
210
211 #endif
212
213 }
214
215 }