Remove unused interactive code from indexutil.*.
authorTeemu Murtola <teemu.murtola@gmail.com>
Wed, 18 Apr 2012 03:46:28 +0000 (06:46 +0300)
committerTeemu Murtola <teemu.murtola@gmail.com>
Wed, 18 Apr 2012 03:46:28 +0000 (06:46 +0300)
The code was a leftover from times when index groups selection was not
handled internally by the selection parser.
Also some minor updates to comments in the header and some formatting
fixes.

Change-Id: Ie5486178011072359ccce884ef7f73fff34a4e55

src/gromacs/selection/indexutil.cpp
src/gromacs/selection/indexutil.h

index dd843eb276c0c04bb13a4d76a72b517bdbf4beb0..c43f150ea790839cf3d3140956ead788e557e9b6 100644 (file)
 #include <config.h>
 #endif
 
-#include <index.h>
-#include <smalloc.h>
-#include <string2.h>
-#include <typedefs.h>
-#include <gmx_fatal.h>
+#include "gromacs/legacyheaders/index.h"
+#include "gromacs/legacyheaders/gmx_fatal.h"
+#include "gromacs/legacyheaders/smalloc.h"
+#include "gromacs/legacyheaders/string2.h"
+#include "gromacs/legacyheaders/typedefs.h"
 
 #include "gromacs/selection/indexutil.h"
 
@@ -74,34 +74,6 @@ gmx_ana_indexgrps_alloc(gmx_ana_indexgrps_t **g, int ngrps)
     snew((*g)->g,    ngrps);
 }
 
-/*!
- * \param[out] g     Index group structure.
- * \param[in]  ngrps Number of index groups.
- * \param[in]  isize Array of index group sizes.
- * \param[in]  index Array of pointers to indices of each group.
- * \param[in]  name  Array of names of the groups.
- * \param[in]  bFree If true, the \p isize, \p index and \p name arrays
- *   are freed after they have been copied.
- */
-void
-gmx_ana_indexgrps_set(gmx_ana_indexgrps_t **g, int ngrps, int *isize,
-                      atom_id **index, char **name, bool bFree)
-{
-    int  i;
-
-    gmx_ana_indexgrps_alloc(g, ngrps);
-    for (i = 0; i < ngrps; ++i)
-    {
-        gmx_ana_index_set(&(*g)->g[i], isize[i], index[i], name[i], isize[i]);
-    }
-    if (bFree)
-    {
-        sfree(isize);
-        sfree(index);
-        sfree(name);
-    }
-}
-
 /*!
  * \param[out] g     Index group structure.
  * \param[in]  top   Topology structure.
@@ -116,7 +88,7 @@ gmx_ana_indexgrps_set(gmx_ana_indexgrps_t **g, int ngrps, int *isize,
  * If both are null, the index group structure is initialized empty.
  */
 void
-gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top, 
+gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
                        const char *fnm)
 {
     t_blocka *block = NULL;
@@ -160,56 +132,6 @@ gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
     sfree(names);
 }
 
-/*!
- * \param[out] g     Index group structure.
- * \param[in]  top   Topology structure.
- * \param[in]  fnm   File name for the index file.
- * \param[in]  ngrps Number of required groups.
- *   Memory is automatically allocated.
- *
- * One of \p top or \p fnm can be NULL, but not both.
- * If \p top is NULL, an index file is required and the groups are read
- * from the file (uses Gromacs routine rd_index()).
- * If \p fnm is NULL, default groups are constructed based on the
- * topology (uses Gromacs routine get_index()).
- */
-void
-gmx_ana_indexgrps_get(gmx_ana_indexgrps_t **g, t_topology *top, 
-                      const char *fnm, int ngrps)
-{
-    int      *isize;
-    atom_id **index;
-    char    **name;
-
-    snew(isize, ngrps);
-    snew(index, ngrps);
-    snew(name,  ngrps);
-    if (!top)
-    {
-        rd_index(fnm, ngrps, isize, index, name);
-    }
-    else
-    {
-        get_index(&(top->atoms), fnm, ngrps, isize, index, name);
-    }
-    gmx_ana_indexgrps_set(g, ngrps, isize, index, name, true);
-}
-
-/*!
- * \param[out] g     Index group structure.
- * \param[in]  fnm   File name for the index file.
- * \param[in]  ngrps Number of required groups.
- *   Memory is automatically allocated.
- *
- * This is a convenience function for calling the Gromacs routine
- * rd_index().
- */
-void
-gmx_ana_indexgrps_rd(gmx_ana_indexgrps_t **g, const char *fnm, int ngrps)
-{
-    gmx_ana_indexgrps_get(g, NULL, fnm, ngrps);
-}
-
 /*!
  * \param[in] g  Index groups structure.
  *
@@ -740,7 +662,6 @@ gmx_ana_index_difference_size(gmx_ana_index_t *a, gmx_ana_index_t *b)
 void
 gmx_ana_index_partition(gmx_ana_index_t *dest1, gmx_ana_index_t *dest2,
                         gmx_ana_index_t *src, gmx_ana_index_t *g)
-                     
 {
     int i, j, k;
 
index 031bfce6cdc74fdb4321aad85478fd88b2d51a5f..1f4e17d703871b4d92338f703bb295a064c88db3 100644 (file)
  * needs to have a unique ID for each possible atom/residue/molecule in the
  * selection, e.g., for analysis of dynamics or for look-up tables.
  *
- * Mostly, these functions are used internally by the library and the
- * selection engine.
- * However, some of the checking functions can be useful in user code to
+ * Mostly, these functions are used internally by the selection engine, but
+ * it is necessary to use some of these functions in order to provide external
+ * index groups to a gmx::SelectionCollection.
+ * Some of the checking functions can be useful outside the selection engine to
  * check the validity of input groups.
- * Also, the mapping functions are useful when dealing with dynamic index
- * groups.
  *
  * \author Teemu Murtola <teemu.murtola@cbr.su.se>
+ * \inlibraryapi
  * \ingroup module_selection
  */
 #ifndef GMX_SELECTION_INDEXUTIL_H
@@ -190,22 +190,10 @@ typedef struct gmx_ana_indexmap_t
 /** Allocate memory for index groups. */
 void
 gmx_ana_indexgrps_alloc(gmx_ana_indexgrps_t **g, int ngrps);
-/** Initializes index groups from arrays. */
-void
-gmx_ana_indexgrps_set(gmx_ana_indexgrps_t **g, int ngrps, int *isize,
-                      atom_id **index, char **name, bool bFree);
 /** Reads index groups from a file or constructs them from topology. */
 void
-gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top, 
+gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
                        const char *fnm);
-/** Ask user to select index groups, possibly constructing groups from 
- *  topology. */
-void
-gmx_ana_indexgrps_get(gmx_ana_indexgrps_t **g, t_topology *top, 
-                      const char *fnm, int ngrps);
-/** Ask user to select index groups from those specified in a file. */
-void
-gmx_ana_indexgrps_rd(gmx_ana_indexgrps_t **g, const char *fnm, int ngrps);
 /** Frees memory allocated for index groups. */
 void
 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t *g);