Move pargs.cpp into src/gromacs/commandline/
authorTeemu Murtola <teemu.murtola@gmail.com>
Sun, 15 Dec 2013 19:43:50 +0000 (21:43 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Sat, 21 Dec 2013 11:53:39 +0000 (12:53 +0100)
- Move the source file and clean it up a bit.
- Move the declarations from readinp.h into a separate (new) header,
  named the same as the source file.

Change-Id: If1775a7387a164bed4fdb2a1d67457d8068f5c6a

src/gromacs/commandline/CMakeLists.txt
src/gromacs/commandline/pargs.cpp [moved from src/gromacs/gmxlib/pargs.cpp with 95% similarity]
src/gromacs/commandline/pargs.h [new file with mode: 0644]
src/gromacs/commandline/wman.h
src/gromacs/legacyheaders/readinp.h
src/gromacs/legacyheaders/statutil.h

index 698c8219e9ff1c63762c510b751d9434df7bcfba..a288ef34b68d8f1d579d7cd1498fb10094b09db9 100644 (file)
@@ -37,7 +37,8 @@ set(LIBGROMACS_SOURCES ${LIBGROMACS_SOURCES} ${COMMANDLINE_SOURCES} PARENT_SCOPE
 
 set(COMMANDLINE_PUBLIC_HEADERS
     cmdlinehelpwriter.h
-    cmdlineparser.h)
+    cmdlineparser.h
+    pargs.h)
 gmx_install_headers(commandline ${COMMANDLINE_PUBLIC_HEADERS})
 
 # mdrun build does not include the selection machinery, which is used in the
similarity index 95%
rename from src/gromacs/gmxlib/pargs.cpp
rename to src/gromacs/commandline/pargs.cpp
index 7ff83218c14849abbe2b209ebdae20f5c0ead1a4..0f965c35af5e3ca6f2a0723f7cd8548d4569a87e 100644 (file)
  * GROningen Mixture of Alchemy and Childrens' Stories
  */
 /* This file is completely threadsafe - keep it that way! */
+#include "gromacs/commandline/pargs.h"
+
 #ifdef HAVE_CONFIG_H
 #include <config.h>
 #endif
 
-#include <ctype.h>
-#include <stdio.h>
-#include <stdlib.h>
-#include <string.h>
-
-#include "typedefs.h"
-#include "gmx_fatal.h"
-#include "statutil.h"
-#include "readinp.h"
-#include "smalloc.h"
-#include "names.h"
-#include "string2.h"
-#include "vec.h"
-#include "macros.h"
+#include <cctype>
+#include <cstdio>
+#include <cstdlib>
+#include <cstring>
+
+#include "gromacs/legacyheaders/gmx_fatal.h"
+#include "gromacs/legacyheaders/smalloc.h"
+#include "gromacs/legacyheaders/string2.h"
 
 #include "gromacs/utility/gmxassert.h"
 
@@ -221,7 +217,7 @@ void get_pargs(int *argc, char *argv[], int nparg, t_pargs pa[], gmx_bool bKeepA
                         *(pa[j].u.c) = sscan(*argc, argv, &i);
                         break;
                     case etENUM:
-                        match = NOTSET;
+                        match = -1;
                         ptr   = sscan(*argc, argv, &i);
                         for (k = 1; (pa[j].u.c[k] != NULL); k++)
                         {
@@ -229,7 +225,7 @@ void get_pargs(int *argc, char *argv[], int nparg, t_pargs pa[], gmx_bool bKeepA
                                pa[j].u.c[k] */
                             if (gmx_strncasecmp(ptr, pa[j].u.c[k], strlen(ptr)) == 0)
                             {
-                                if ( ( match == NOTSET ) ||
+                                if ( ( match == -1 ) ||
                                      ( strlen(pa[j].u.c[k]) <
                                        strlen(pa[j].u.c[match]) ) )
                                 {
@@ -237,7 +233,7 @@ void get_pargs(int *argc, char *argv[], int nparg, t_pargs pa[], gmx_bool bKeepA
                                 }
                             }
                         }
-                        if (match != NOTSET)
+                        if (match != -1)
                         {
                             pa[j].u.c[0] = pa[j].u.c[match];
                         }
@@ -325,7 +321,7 @@ gmx_bool opt2parg_gmx_bool(const char *option, int nparg, t_pargs pa[])
         }
     }
 
-    gmx_fatal(FARGS, "No gmx_boolean option %s in pargs", option);
+    gmx_fatal(FARGS, "No boolean option %s in pargs", option);
 
     return FALSE;
 }
diff --git a/src/gromacs/commandline/pargs.h b/src/gromacs/commandline/pargs.h
new file mode 100644 (file)
index 0000000..04e17cd
--- /dev/null
@@ -0,0 +1,105 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
+ * Copyright (c) 2001-2004, The GROMACS development team.
+ * Copyright (c) 2013, 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.
+ */
+#ifndef GMX_COMMANDLINE_PARGS_H
+#define GMX_COMMANDLINE_PARGS_H
+
+#include "../legacyheaders/types/simple.h"
+
+#ifdef __cplusplus
+extern "C"
+{
+#endif
+
+/* This structure is used for parsing arguments off the comand line */
+enum
+{
+    etINT, etGMX_LARGE_INT, etREAL, etTIME, etSTR,    etBOOL, etRVEC,   etENUM, etNR
+};
+
+typedef struct
+{
+    const char *option;
+    gmx_bool    bSet;
+    int         type;
+    union
+    {
+        void            *v; /* This is a nasty workaround, to be able to use initialized */
+        int             *i; /* arrays */
+        gmx_large_int_t *is;
+        real            *r;
+        const char     **c; /* Must be pointer to string (when type == etSTR)         */
+        /* or null terminated list of enums (when type == etENUM) */
+        gmx_bool        *b;
+        rvec            *rv;
+    }           u;
+    const char *desc;
+} t_pargs;
+
+void get_pargs(int *argc, char *argv[], int nparg, t_pargs pa[],
+               gmx_bool bKeepArgs);
+/* Read a number of arguments from the command line.
+ * For etINT, etREAL and etCHAR an extra argument is read (when present)
+ * for etBOOL the gmx_boolean option is changed to the negate value
+ * If !bKeepArgs, the command line arguments are removed from the command line
+ */
+
+gmx_bool is_hidden(t_pargs *pa);
+/* Return TRUE when the option is a secret one */
+
+int nenum(const char *const enumc[]);
+/* returns ordinal number of selected enum from args
+ * depends on enumc[0] pointing to one of the other elements
+ * array must be terminated by a NULL pointer
+ */
+
+int opt2parg_int(const char *option, int nparg, t_pargs pa[]);
+
+gmx_bool opt2parg_gmx_bool(const char *option, int nparg, t_pargs pa[]);
+
+real opt2parg_real(const char *option, int nparg, t_pargs pa[]);
+
+const char *opt2parg_str(const char *option, int nparg, t_pargs pa[]);
+
+const char *opt2parg_enum(const char *option, int nparg, t_pargs pa[]);
+
+gmx_bool opt2parg_bSet(const char *option, int nparg, t_pargs pa[]);
+
+#ifdef __cplusplus
+}
+#endif
+
+#endif
index 8e6ba07980538a2e85ef7a2724d8dd8f6c828aed..ec5631451c8b5e71f2193aea622cd6e5fc279090 100644 (file)
@@ -37,8 +37,7 @@
 #ifndef GMX_COMMANDLINE_WMAN_H
 #define GMX_COMMANDLINE_WMAN_H
 
-#include "gromacs/legacyheaders/readinp.h"
-
+#include "gromacs/commandline/pargs.h"
 #include "gromacs/fileio/filenm.h"
 
 void write_man(const char *mantp, const char *program,
index 418ecb0ac36baaf1fc5621020034d67fa6eb0226..31fd5c31cfb02336d651d68f4802a574e2228006 100644 (file)
@@ -39,7 +39,6 @@
 #include "typedefs.h"
 #include "warninp.h"
 
-
 #ifdef __cplusplus
 extern "C" {
 #endif
@@ -57,8 +56,6 @@ typedef struct {
    Initally read in with read_inpfile, then filled in with missing values
    through get_eint, get_ereal, etc. */
 
-
-
 t_inpfile *read_inpfile(const char *fn, int *ninp,
                         warninp_t wi);
 /* Create & populate a t_inpfile struct from values in file fn.
@@ -109,57 +106,6 @@ int get_eenum(int *ninp, t_inpfile **inp, const char *name, const char **defs);
 #define CTYPE(s)  STYPENC("; " s, NULL)
 /* This last one prints a comment line where you can add some explanation */
 
-/* This structure is used for parsing arguments off the comand line */
-enum {
-    etINT, etGMX_LARGE_INT, etREAL, etTIME, etSTR,    etBOOL, etRVEC,   etENUM, etNR
-};
-
-typedef struct {
-    const char *option;
-    gmx_bool    bSet;
-    int         type;
-    union {
-        void            *v; /* This is a nasty workaround, to be able to use initialized */
-        int             *i; /* arrays */
-        gmx_large_int_t *is;
-        real            *r;
-        const char     **c; /* Must be pointer to string (when type == etSTR)         */
-        /* or null terminated list of enums (when type == etENUM) */
-        gmx_bool        *b;
-        rvec            *rv;
-    } u;
-    const char *desc;
-} t_pargs;
-
-void get_pargs(int *argc, char *argv[], int nparg, t_pargs pa[],
-               gmx_bool bKeepArgs);
-/* Read a number of arguments from the command line.
- * For etINT, etREAL and etCHAR an extra argument is read (when present)
- * for etBOOL the gmx_boolean option is changed to the negate value
- * If !bKeepArgs, the command line arguments are removed from the command line
- */
-
-gmx_bool is_hidden(t_pargs *pa);
-/* Return TRUE when the option is a secret one */
-
-int nenum(const char *const enumc[]);
-/* returns ordinal number of selected enum from args
- * depends on enumc[0] pointing to one of the other elements
- * array must be terminated by a NULL pointer
- */
-
-int opt2parg_int(const char *option, int nparg, t_pargs pa[]);
-
-gmx_bool opt2parg_gmx_bool(const char *option, int nparg, t_pargs pa[]);
-
-real opt2parg_real(const char *option, int nparg, t_pargs pa[]);
-
-const char *opt2parg_str(const char *option, int nparg, t_pargs pa[]);
-
-const char *opt2parg_enum(const char *option, int nparg, t_pargs pa[]);
-
-gmx_bool opt2parg_bSet(const char *option, int nparg, t_pargs pa[]);
-
 #ifdef __cplusplus
 }
 #endif
index b6af4be90b847d8861432a550fda045997ba5826..b3f0c7814ea6b915d14e9135bc928c0302a1d638 100644 (file)
@@ -36,8 +36,8 @@
 #ifndef _statutil_h
 #define _statutil_h
 
+#include "../commandline/pargs.h"
 #include "../fileio/filenm.h"
-#include "readinp.h"
 #include "oenv.h"
 
 #ifdef __cplusplus