Skip to content
Snippets Groups Projects
Commit d0071abf authored by Mark Abraham's avatar Mark Abraham Committed by Gerrit Code Review
Browse files

Fix grompp not warning when .mdp values have wrong types

With old code, e.g. init_lambda_state = 0.35 was silently
interpreted as 0.

Added some comments and TODOs.

Added unit tests for some simple parts of the the functionality that
needed to be fixed. Extended the interface of warninp_t to permit
sound testing.

Fixes #1893

Change-Id: I425798daa8ed3a391cde49ba4a2d3d352d682f8f
parent 2ee7caef
No related branches found
No related tags found
Loading
#
# This file is part of the GROMACS molecular simulation package.
#
# Copyright (c) 2009,2010,2012,2014,2015, by the GROMACS development team, led by
# Copyright (c) 2009,2010,2012,2014,2015,2016, by the GROMACS development team, led by
# Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
# and including many others, as listed in the AUTHORS file in the
# top-level source directory and at http://www.gromacs.org.
......@@ -51,3 +51,7 @@ endif()
add_subdirectory(gpu_utils)
set(GMXLIB_SOURCES ${GMXLIB_SOURCES} ${NONBONDED_SOURCES} PARENT_SCOPE)
if(BUILD_TESTING)
add_subdirectory(tests)
endif()
......@@ -3,7 +3,7 @@
*
* Copyright (c) 1991-2000, University of Groningen, The Netherlands.
* Copyright (c) 2001-2004, The GROMACS development team.
* Copyright (c) 2013,2014,2015, by the GROMACS development team, led by
* Copyright (c) 2013,2014,2015,2016, by the GROMACS development team, led by
* Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
* and including many others, as listed in the AUTHORS file in the
* top-level source directory and at http://www.gromacs.org.
......@@ -79,6 +79,9 @@ t_inpfile *read_inpfile(const char *fn, int *ninp,
set_warning_line(wi, fn, lc);
if (ptr)
{
// TODO This parsing should be using strip_comment, trim,
// strchr, etc. rather than re-inventing wheels.
/* Strip comment */
if ((cptr = strchr(buf, COMMENTSIGN)) != NULL)
{
......@@ -366,6 +369,7 @@ static int get_einp(int *ninp, t_inpfile **inp, const char *name)
}
}
/* Note that sanitizing the trailing part of (*inp)[ii].value was the responsibility of read_inpfile() */
int get_eint(int *ninp, t_inpfile **inp, const char *name, int def,
warninp_t wi)
{
......@@ -385,7 +389,7 @@ int get_eint(int *ninp, t_inpfile **inp, const char *name, int def,
else
{
ret = strtol((*inp)[ii].value, &ptr, 10);
if (ptr == (*inp)[ii].value)
if (*ptr != '\0')
{
sprintf(warn_buf, "Right hand side '%s' for parameter '%s' in parameter file is not an integer value\n", (*inp)[ii].value, (*inp)[ii].name);
warning_error(wi, warn_buf);
......@@ -395,6 +399,7 @@ int get_eint(int *ninp, t_inpfile **inp, const char *name, int def,
}
}
/* Note that sanitizing the trailing part of (*inp)[ii].value was the responsibility of read_inpfile() */
gmx_int64_t get_eint64(int *ninp, t_inpfile **inp,
const char *name, gmx_int64_t def,
warninp_t wi)
......@@ -415,7 +420,7 @@ gmx_int64_t get_eint64(int *ninp, t_inpfile **inp,
else
{
ret = str_to_int64_t((*inp)[ii].value, &ptr);
if (ptr == (*inp)[ii].value)
if (*ptr != '\0')
{
sprintf(warn_buf, "Right hand side '%s' for parameter '%s' in parameter file is not an integer value\n", (*inp)[ii].value, (*inp)[ii].name);
warning_error(wi, warn_buf);
......@@ -425,6 +430,7 @@ gmx_int64_t get_eint64(int *ninp, t_inpfile **inp,
}
}
/* Note that sanitizing the trailing part of (*inp)[ii].value was the responsibility of read_inpfile() */
double get_ereal(int *ninp, t_inpfile **inp, const char *name, double def,
warninp_t wi)
{
......@@ -444,7 +450,7 @@ double get_ereal(int *ninp, t_inpfile **inp, const char *name, double def,
else
{
ret = strtod((*inp)[ii].value, &ptr);
if (ptr == (*inp)[ii].value)
if (*ptr != '\0')
{
sprintf(warn_buf, "Right hand side '%s' for parameter '%s' in parameter file is not a real value\n", (*inp)[ii].value, (*inp)[ii].name);
warning_error(wi, warn_buf);
......@@ -454,6 +460,7 @@ double get_ereal(int *ninp, t_inpfile **inp, const char *name, double def,
}
}
/* Note that sanitizing the trailing part of (*inp)[ii].value was the responsibility of read_inpfile() */
const char *get_estr(int *ninp, t_inpfile **inp, const char *name, const char *def)
{
char buf[32];
......@@ -481,6 +488,7 @@ const char *get_estr(int *ninp, t_inpfile **inp, const char *name, const char *d
}
}
/* Note that sanitizing the trailing part of (*inp)[ii].value was the responsibility of read_inpfile() */
int get_eeenum(int *ninp, t_inpfile **inp, const char *name, const char **defs,
warninp_t wi)
{
......
#
# This file is part of the GROMACS molecular simulation package.
#
# Copyright (c) 2016, by the GROMACS development team, led by
# Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
# and including many others, as listed in the AUTHORS file in the
# top-level source directory and at http://www.gromacs.org.
#
# GROMACS is free software; you can redistribute it and/or
# modify it under the terms of the GNU Lesser General Public License
# as published by the Free Software Foundation; either version 2.1
# of the License, or (at your option) any later version.
#
# GROMACS is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
# Lesser General Public License for more details.
#
# You should have received a copy of the GNU Lesser General Public
# License along with GROMACS; if not, see
# http://www.gnu.org/licenses, or write to the Free Software Foundation,
# Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
#
# If you want to redistribute modifications to GROMACS, please
# consider that scientific software is very special. Version
# control is crucial - bugs must be traceable. We will be happy to
# consider code for inclusion in the official distribution, but
# derived work must not be called official GROMACS. Details are found
# in the README & COPYING files - if they are missing, get the
# official version at http://www.gromacs.org.
#
# To help us fund GROMACS development, we humbly ask that you cite
# the research papers on the package. Check out http://www.gromacs.org.
gmx_add_unit_test(GmxlibTests gmxlib-test
readinp.cpp
)
/*
* This file is part of the GROMACS molecular simulation package.
*
* Copyright (c) 2016, by the GROMACS development team, led by
* Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
* and including many others, as listed in the AUTHORS file in the
* top-level source directory and at http://www.gromacs.org.
*
* GROMACS is free software; you can redistribute it and/or
* modify it under the terms of the GNU Lesser General Public License
* as published by the Free Software Foundation; either version 2.1
* of the License, or (at your option) any later version.
*
* GROMACS is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
* Lesser General Public License for more details.
*
* You should have received a copy of the GNU Lesser General Public
* License along with GROMACS; if not, see
* http://www.gnu.org/licenses, or write to the Free Software Foundation,
* Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
*
* If you want to redistribute modifications to GROMACS, please
* consider that scientific software is very special. Version
* control is crucial - bugs must be traceable. We will be happy to
* consider code for inclusion in the official distribution, but
* derived work must not be called official GROMACS. Details are found
* in the README & COPYING files - if they are missing, get the
* official version at http://www.gromacs.org.
*
* To help us fund GROMACS development, we humbly ask that you cite
* the research papers on the package. Check out http://www.gromacs.org.
*/
/*! \internal \file
* \brief
* Tests utilities for routines that parse fields e.g. from grompp input
*
* \author Mark Abraham <mark.j.abraham@gmail.com>
*/
#include "gmxpre.h"
#include "gromacs/legacyheaders/readinp.h"
#include <gtest/gtest.h>
#include "gromacs/legacyheaders/warninp.h"
#include "gromacs/utility/scoped_cptr.h"
#include "gromacs/utility/smalloc.h"
namespace gmx
{
namespace testing
{
class ReadTest : public ::testing::Test
{
public:
ReadTest() : numInputs_(1),
inputField_(0),
inpGuard_(),
wi_(),
wiGuard_()
{
snew(inputField_, numInputs_);
inpGuard_.reset(inputField_);
inputField_[0].count = 0;
inputField_[0].bObsolete = FALSE;
inputField_[0].bSet = FALSE;
inputField_[0].name = (char *) "test";
inputField_[0].inp_count = 0;
wi_ = init_warning(FALSE, 0);
wiGuard_.reset(wi_);
}
int numInputs_;
t_inpfile *inputField_;
gmx::scoped_cptr<t_inpfile> inpGuard_;
warninp_t wi_;
gmx::scoped_cptr<struct warninp, free_warning> wiGuard_;
};
TEST_F(ReadTest, get_eint_ReadsInteger)
{
inputField_[0].value = (char *) "1";
ASSERT_EQ(1, get_eint(&numInputs_, &inputField_, "test", 2, wi_));
ASSERT_FALSE(warning_errors_exist(wi_));
}
TEST_F(ReadTest, get_eint_WarnsAboutFloat)
{
inputField_[0].value = (char *) "0.8";
get_eint(&numInputs_, &inputField_, "test", 2, wi_);
ASSERT_TRUE(warning_errors_exist(wi_));
}
TEST_F(ReadTest, get_eint_WarnsAboutString)
{
inputField_[0].value = (char *) "hello";
get_eint(&numInputs_, &inputField_, "test", 2, wi_);
ASSERT_TRUE(warning_errors_exist(wi_));
}
TEST_F(ReadTest, get_eint64_ReadsInteger)
{
inputField_[0].value = (char *) "1";
ASSERT_EQ(1, get_eint64(&numInputs_, &inputField_, "test", 2, wi_));
ASSERT_FALSE(warning_errors_exist(wi_));
}
TEST_F(ReadTest, get_eint64_WarnsAboutFloat)
{
inputField_[0].value = (char *) "0.8";
get_eint64(&numInputs_, &inputField_, "test", 2, wi_);
ASSERT_TRUE(warning_errors_exist(wi_));
}
TEST_F(ReadTest, get_eint64_WarnsAboutString)
{
inputField_[0].value = (char *) "hello";
get_eint64(&numInputs_, &inputField_, "test", 2, wi_);
ASSERT_TRUE(warning_errors_exist(wi_));
}
TEST_F(ReadTest, get_ereal_ReadsInteger)
{
inputField_[0].value = (char *) "1";
ASSERT_EQ(1, get_ereal(&numInputs_, &inputField_, "test", 2, wi_));
ASSERT_FALSE(warning_errors_exist(wi_));
}
TEST_F(ReadTest, get_ereal_ReadsFloat)
{
inputField_[0].value = (char *) "0.8";
ASSERT_EQ(0.8, get_ereal(&numInputs_, &inputField_, "test", 2, wi_));
ASSERT_FALSE(warning_errors_exist(wi_));
}
TEST_F(ReadTest, get_ereal_WarnsAboutString)
{
inputField_[0].value = (char *) "hello";
get_ereal(&numInputs_, &inputField_, "test", 2, wi_);
ASSERT_TRUE(warning_errors_exist(wi_));
}
} // namespace
} // namespace
......@@ -3,7 +3,7 @@
*
* Copyright (c) 1991-2000, University of Groningen, The Netherlands.
* Copyright (c) 2001-2004, The GROMACS development team.
* Copyright (c) 2013,2014,2015, by the GROMACS development team, led by
* Copyright (c) 2013,2014,2015,2016, by the GROMACS development team, led by
* Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
* and including many others, as listed in the AUTHORS file in the
* top-level source directory and at http://www.gromacs.org.
......@@ -177,6 +177,11 @@ void check_warning_error(warninp_t wi, int f_errno, const char *file, int line)
}
}
gmx_bool warning_errors_exist(warninp_t wi)
{
return (wi->nwarn_error > 0);
}
void done_warning(warninp_t wi, int f_errno, const char *file, int line)
{
print_warn_count("note", wi->nwarn_note);
......@@ -192,6 +197,11 @@ void done_warning(warninp_t wi, int f_errno, const char *file, int line)
wi->nwarn_warn, ShortProgram());
}
free_warning(wi);
}
void free_warning(warninp_t wi)
{
sfree(wi);
}
......
......@@ -3,7 +3,7 @@
*
* Copyright (c) 1991-2000, University of Groningen, The Netherlands.
* Copyright (c) 2001-2004, The GROMACS development team.
* Copyright (c) 2010,2014, by the GROMACS development team, led by
* Copyright (c) 2010,2014,2016, by the GROMACS development team, led by
* Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
* and including many others, as listed in the AUTHORS file in the
* top-level source directory and at http://www.gromacs.org.
......@@ -101,6 +101,9 @@ warning_error(warninp_t wi, const char *s);
* are printed, nwarn_error (local) is incremented.
*/
gmx_bool warning_errors_exist(warninp_t wi);
/* Return whether any error-level warnings were issued to wi. */
void
check_warning_error(warninp_t wi, int f_errno, const char *file, int line);
/* When warning_error has been called at least once gmx_fatal is called,
......@@ -116,6 +119,10 @@ done_warning(warninp_t wi, int f_errno, const char *file, int line);
* Frees the data structure pointed to by wi.
*/
void
free_warning(warninp_t wi);
/* Frees the data structure pointed to by wi. */
void
_too_few(warninp_t wi, const char *fn, int line);
#define too_few(wi) _too_few(wi, __FILE__, __LINE__)
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment