Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / mdlib / groupcoord.h
index 4893ba130bcbd65c2650ac18c9ee72e980191a35..686ce17517e2a1cc82c89801570fcd4383fc1158 100644 (file)
  * To help us fund GROMACS development, we humbly ask that you cite
  * the research papers on the package. Check out http://www.gromacs.org.
  */
-/*! \libinternal \file groupcoord.h
- * \brief
- * Assemble atom positions for comparison with a reference set.
+
+/*! \libinternal \file
+ * \brief Assemble atomic positions of a (small) subset of atoms and distribute to all nodes.
  *
- * This file contains functions to assemble the positions of a subset of the
- * atoms and to do operations on it like determining the center of mass, or
- * doing translations and rotations. These functions are useful when
- * a subset of the positions needs to be compared to some set of reference
- * positions, as e.g. done for essential dynamics.
+ *  This file contains functions to assemble the positions of a subset of the
+ *  atoms and to do operations on it like determining the center of mass, or
+ *  doing translations and rotations. These functions are useful when
+ *  a subset of the positions needs to be compared to some set of reference
+ *  positions, as e.g. done for essential dynamics.
  *
  * \inlibraryapi
  */
 
-#ifdef HAVE_CONFIG_H
-#include <config.h>
-#endif
-
 #include <stdio.h>
-#include "typedefs.h"
-#include "types/commrec.h"
 
+#include "gromacs/legacyheaders/typedefs.h"
+#include "gromacs/legacyheaders/types/commrec.h"
+
+#ifdef __cplusplus
+extern "C" {
+#endif
 
 /*! \brief Select local atoms of a group.
  *
@@ -74,6 +74,7 @@
  *                            (collective) array such that it can be gmx_summed
  *                            in the communicate_group_positions routine.
  */
+
 extern void dd_make_local_group_indices(gmx_ga2la_t ga2la,
                                         const int nr, int anrs[], int *nr_loc,
                                         int *anrs_loc[], int *nalloc_loc,
@@ -88,19 +89,28 @@ extern void dd_make_local_group_indices(gmx_ga2la_t ga2la,
  * retrieved from anrs_loc[0..nr_loc]. If you call the routine for the serial case,
  * provide an array coll_ind[i] = i for i in 1..nr.
  *
+ * If shifts != NULL, the PBC representation of each atom is chosen such that a
+ * continuous trajectory results. Therefore, if the group is whole at the start
+ * of the simulation, it will always stay whole.
+ * If shifts = NULL, the group positions are not made whole again, but assembled
+ * and distributed to all nodes. The variables marked "optional" are not used in
+ * that case.
+ *
  * \param[in]     cr           Pointer to MPI communication data.
- * \param[out]    xcoll        Collective array of positions, idential on all nodes
+ * \param[out]    xcoll        Collective array of positions, identical on all nodes
  *                             after this routine has been called.
  * \param[in,out] shifts       Collective array of shifts for xcoll, needed to make
  *                             the group whole. This array remembers the shifts
  *                             since the start of the simulation (where the group
  *                             is whole) and must therefore not be changed outside
- *                             of this routine!
+ *                             of this routine! If NULL, the group will not be made
+ *                             whole and the optional variables are ignored.
  * \param[out]    extra_shifts Extra shifts since last time step, only needed as
- *                             buffer variable [0..nr].
- * \param[in]     bNS          Neighborsearching/domain redecomposition has been
+ *                             buffer variable [0..nr] (optional).
+ * \param[in]     bNS          Neighbor searching / domain re-decomposition has been
  *                             performed at the begin of this time step such that
- *                             the shifts have changed and need to be updated.
+ *                             the shifts have changed and need to be updated
+ *                             (optional).
  * \param[in]     x_loc        Pointer to the local atom positions this node has.
  * \param[in]     nr           Total number of atoms in the group.
  * \param[in]     nr_loc       Number of group atoms on the local node.
@@ -110,9 +120,9 @@ extern void dd_make_local_group_indices(gmx_ga2la_t ga2la,
  *                             contributions can be gmx_summed. It is provided by
  *                             dd_make_local_group_indices.
  * \param[in,out] xcoll_old    Positions from the last time step, used to make the
- *                             group whole.
+ *                             group whole (optional).
  * \param[in]     box          Simulation box matrix, needed to shift xcoll such that
- *                             the group becomes whole.
+ *                             the group becomes whole (optional).
  */
 extern void communicate_group_positions(t_commrec *cr, rvec *xcoll, ivec *shifts,
                                         ivec *extra_shifts, const gmx_bool bNS,
@@ -120,7 +130,6 @@ extern void communicate_group_positions(t_commrec *cr, rvec *xcoll, ivec *shifts
                                         int *anrs_loc, int *coll_ind, rvec *xcoll_old,
                                         matrix box);
 
-
 /*! \brief Calculates the center of the positions x locally.
  *
  * Calculates the center of mass (if masses are given in the weight array) or
@@ -196,3 +205,7 @@ extern void translate_x(rvec x[], const int nr, const rvec transvec);
  *
  */
 extern void rotate_x(rvec x[], const int nr, matrix rmat);
+
+#ifdef __cplusplus
+}
+#endif