Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / selection / sm_permute.cpp
index e6b86cf543e3cae65599403b53c9db25348fbb42..a4f52fe150e3e97b472a5e59996b4446153778f5 100644 (file)
@@ -1,49 +1,55 @@
 /*
+ * This file is part of the GROMACS molecular simulation package.
  *
- *                This source code is part of
+ * Copyright (c) 2009,2010,2011,2012,2013,2014, 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.
  *
- *                 G   R   O   M   A   C   S
+ * 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.
  *
- *          GROningen MAchine for Chemical Simulations
+ * 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.
  *
- * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
- * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
- * Copyright (c) 2001-2009, The GROMACS development team,
- * check out http://www.gromacs.org for more information.
-
- * This program is free software; you can redistribute it and/or
- * modify it under the terms of the GNU General Public License
- * as published by the Free Software Foundation; either version 2
- * of the License, or (at your option) any later version.
+ * 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, 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 www.gromacs.org.
+ * 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 papers on the package - you can find them in the top README file.
- *
- * For more info, check our website at http://www.gromacs.org
+ * the research papers on the package. Check out http://www.gromacs.org.
  */
 /*! \internal \file
  * \brief
  * Implements the \p permute selection modifier.
  *
- * \author Teemu Murtola <teemu.murtola@cbr.su.se>
+ * \author Teemu Murtola <teemu.murtola@gmail.com>
  * \ingroup module_selection
  */
-#include "gromacs/legacyheaders/macros.h"
-#include "gromacs/legacyheaders/smalloc.h"
-#include "gromacs/legacyheaders/vec.h"
+#include "gmxpre.h"
 
+#include "gromacs/legacyheaders/macros.h"
+#include "gromacs/math/vec.h"
 #include "gromacs/selection/position.h"
-#include "gromacs/selection/selmethod.h"
 #include "gromacs/utility/exceptions.h"
+#include "gromacs/utility/smalloc.h"
 #include "gromacs/utility/stringutil.h"
 
+#include "selmethod.h"
+
 /*! \internal \brief
  * Data structure for the \p permute selection modifier.
  */
@@ -51,8 +57,6 @@ typedef struct
 {
     /** Positions to permute. */
     gmx_ana_pos_t    p;
-    /** Group to receive the output permutation. */
-    gmx_ana_index_t  g;
     /** Number of elements in the permutation. */
     int              n;
     /** Array describing the permutation. */
@@ -61,20 +65,56 @@ typedef struct
     int             *rperm;
 } t_methoddata_permute;
 
-/** Allocates data for the \p permute selection modifier. */
+/*! \brief
+ * Allocates data for the \p permute selection modifier.
+ *
+ * \param[in]     npar  Not used (should be 2).
+ * \param[in,out] param Method parameters (should point to a copy of
+ *   \ref smparams_permute).
+ * \returns Pointer to the allocated data (\p t_methoddata_permute).
+ *
+ * Allocates memory for a \p t_methoddata_permute structure.
+ */
 static void *
 init_data_permute(int npar, gmx_ana_selparam_t *param);
-/** Initializes data for the \p permute selection modifier. */
+/*! \brief
+ * Initializes data for the \p permute selection modifier.
+ *
+ * \param[in] top   Not used.
+ * \param[in] npar  Not used (should be 2).
+ * \param[in] param Method parameters (should point to \ref smparams_permute).
+ * \param[in] data  Should point to a \p t_methoddata_permute.
+ * \returns   0 if the input permutation is valid, -1 on error.
+ */
 static void
 init_permute(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
-/** Initializes output for the \p permute selection modifier. */
+/*! \brief
+ * Initializes output for the \p permute selection modifier.
+ *
+ * \param[in]     top   Topology data structure.
+ * \param[in,out] out   Pointer to output data structure.
+ * \param[in,out] data  Should point to \c t_methoddata_permute.
+ */
 static void
 init_output_permute(t_topology *top, gmx_ana_selvalue_t *out, void *data);
 /** Frees the memory allocated for the \p permute selection modifier. */
 static void
 free_data_permute(void *data);
-/** Evaluates the \p permute selection modifier. */
 static void
+/*! \brief
+ * Evaluates the \p permute selection modifier.
+ *
+ * \param[in]  top   Not used.
+ * \param[in]  fr    Not used.
+ * \param[in]  pbc   Not used.
+ * \param[in]  p     Positions to permute (should point to \p data->p).
+ * \param[out] out   Output data structure (\p out->u.p is used).
+ * \param[in]  data  Should point to a \p t_methoddata_permute.
+ * \returns    0 if \p p could be permuted, -1 on error.
+ *
+ * Returns -1 if the size of \p p is not divisible by the number of
+ * elements in the permutation.
+ */
 evaluate_permute(t_topology *top, t_trxframe *fr, t_pbc *pbc,
                  gmx_ana_pos_t *p, gmx_ana_selvalue_t *out, void *data);
 
@@ -103,7 +143,7 @@ static const char *help_permute[] = {
     "a selection.",
 };
 
-/** \internal Selection method data for the \p permute modifier. */
+/** Selection method data for the \p permute modifier. */
 gmx_ana_selmethod_t sm_permute = {
     "permute", POS_VALUE, SMETH_MODIFIER,
     asize(smparams_permute), smparams_permute,
@@ -118,44 +158,29 @@ gmx_ana_selmethod_t sm_permute = {
     {"permute P1 ... PN", asize(help_permute), help_permute},
 };
 
-/*!
- * \param[in]     npar  Not used (should be 2).
- * \param[in,out] param Method parameters (should point to a copy of
- *   \ref smparams_permute).
- * \returns Pointer to the allocated data (\p t_methoddata_permute).
- *
- * Allocates memory for a \p t_methoddata_permute structure.
- */
 static void *
-init_data_permute(int npar, gmx_ana_selparam_t *param)
+init_data_permute(int /* npar */, gmx_ana_selparam_t *param)
 {
-    t_methoddata_permute *data;
-
-    snew(data, 1);
+    t_methoddata_permute *data = new t_methoddata_permute();
+    data->n          = 0;
+    data->perm       = NULL;
+    data->rperm      = NULL;
     param[0].val.u.p = &data->p;
     return data;
 }
 
-/*!
- * \param[in] top   Not used.
- * \param[in] npar  Not used (should be 2).
- * \param[in] param Method parameters (should point to \ref smparams_permute).
- * \param[in] data  Should point to a \p t_methoddata_permute.
- * \returns   0 if the input permutation is valid, -1 on error.
- */
 static void
-init_permute(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
+init_permute(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t *param, void *data)
 {
     t_methoddata_permute *d = (t_methoddata_permute *)data;
     int                   i;
 
-    gmx_ana_index_reserve(&d->g, d->p.g->isize);
     d->n    = param[1].val.nr;
     d->perm = param[1].val.u.i;
-    if (d->p.nr % d->n != 0)
+    if (d->p.count() % d->n != 0)
     {
         GMX_THROW(gmx::InconsistentInputError(
-                    gmx::formatString("The number of positions to be permuted is not divisible by %d", d->n)));
+                          gmx::formatString("The number of positions to be permuted is not divisible by %d", d->n)));
     }
     snew(d->rperm, d->n);
     for (i = 0; i < d->n; ++i)
@@ -177,27 +202,22 @@ init_permute(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
     }
 }
 
-/*!
- * \param[in]     top   Topology data structure.
- * \param[in,out] out   Pointer to output data structure.
- * \param[in,out] data  Should point to \c t_methoddata_permute.
- */
 static void
-init_output_permute(t_topology *top, gmx_ana_selvalue_t *out, void *data)
+init_output_permute(t_topology * /* top */, gmx_ana_selvalue_t *out, void *data)
 {
     t_methoddata_permute *d = (t_methoddata_permute *)data;
     int                   i, j, b;
 
-    gmx_ana_pos_copy(out->u.p, &d->p, true);
-    gmx_ana_pos_set_evalgrp(out->u.p, &d->g);
-    d->g.isize = 0;
+    out->u.p->m.type = d->p.m.type;
+    gmx_ana_pos_reserve_for_append(out->u.p, d->p.count(), d->p.m.b.nra,
+                                   d->p.v != NULL, d->p.f != NULL);
     gmx_ana_pos_empty_init(out->u.p);
-    for (i = 0; i < d->p.nr; i += d->n)
+    for (i = 0; i < d->p.count(); i += d->n)
     {
         for (j = 0; j < d->n; ++j)
         {
             b = i + d->rperm[j];
-            gmx_ana_pos_append_init(out->u.p, &d->g, &d->p, b);
+            gmx_ana_pos_append_init(out->u.p, &d->p, b);
         }
     }
 }
@@ -212,50 +232,36 @@ free_data_permute(void *data)
 {
     t_methoddata_permute *d = (t_methoddata_permute *)data;
 
-    gmx_ana_index_deinit(&d->g);
     sfree(d->rperm);
-    sfree(d);
+    delete d;
 }
 
-/*!
- * \param[in]  top   Not used.
- * \param[in]  fr    Not used.
- * \param[in]  pbc   Not used.
- * \param[in]  p     Positions to permute (should point to \p data->p).
- * \param[out] out   Output data structure (\p out->u.p is used).
- * \param[in]  data  Should point to a \p t_methoddata_permute.
- * \returns    0 if \p p could be permuted, -1 on error.
- *
- * Returns -1 if the size of \p p is not divisible by the number of
- * elements in the permutation.
- */
 static void
-evaluate_permute(t_topology *top, t_trxframe *fr, t_pbc *pbc,
+evaluate_permute(t_topology * /* top */, t_trxframe * /* fr */, t_pbc * /* pbc */,
                  gmx_ana_pos_t *p, gmx_ana_selvalue_t *out, void *data)
 {
     t_methoddata_permute *d = (t_methoddata_permute *)data;
     int                   i, j, b;
     int                   refid;
 
-    if (d->p.nr % d->n != 0)
+    if (d->p.count() % d->n != 0)
     {
         GMX_THROW(gmx::InconsistentInputError(
-                    gmx::formatString("The number of positions to be permuted is not divisible by %d", d->n)));
+                          gmx::formatString("The number of positions to be permuted is not divisible by %d", d->n)));
     }
-    d->g.isize = 0;
     gmx_ana_pos_empty(out->u.p);
-    for (i = 0; i < d->p.nr; i += d->n)
+    for (i = 0; i < d->p.count(); i += d->n)
     {
         for (j = 0; j < d->n; ++j)
         {
-            b = i + d->rperm[j];
+            b     = i + d->rperm[j];
             refid = d->p.m.refid[b];
             if (refid != -1)
             {
                 /* De-permute the reference ID */
                 refid = refid - (refid % d->n) + d->perm[refid % d->n];
             }
-            gmx_ana_pos_append(out->u.p, &d->g, p, b, refid);
+            gmx_ana_pos_append(out->u.p, p, b, refid);
         }
     }
     gmx_ana_pos_append_finish(out->u.p);