diff --git a/.gitignore b/.gitignore
index 8c3ede8f3125f5024c1fe01a405024d1cf5a7f19..c31cac8ed4efffdc75fa91f4b608265379c6d48c 100644
--- a/.gitignore
+++ b/.gitignore
@@ -80,6 +80,7 @@ tests/test_nonsym_force_1_vec.dat
 tests/test_nonsym_force_2_vec.dat
 tests/potential.dat
 tests/testGreetings
+tests/testSelectOutput
 tests/testReading
 tests/testSingle
 tests/testTimeIntegration
@@ -102,6 +103,7 @@ tests/test125cells.sh
 tests/test125cellsPerturbed.sh
 tests/testParser.sh
 tests/testReading.sh
+tests/testSelectOutput.sh
 tests/testAdiabaticIndex
 tests/testRiemannExact
 tests/testRiemannTRRS
diff --git a/configure.ac b/configure.ac
index 316155b4962cb853d4c3bb51e5ffbf6c1c0e038e..48a36893b8926224679034ddca0f1ebccc7d275d 100644
--- a/configure.ac
+++ b/configure.ac
@@ -1275,6 +1275,7 @@ AC_CONFIG_FILES([tests/testPeriodicBC.sh], [chmod +x tests/testPeriodicBC.sh])
 AC_CONFIG_FILES([tests/testPeriodicBCPerturbed.sh], [chmod +x tests/testPeriodicBCPerturbed.sh])
 AC_CONFIG_FILES([tests/testInteractions.sh], [chmod +x tests/testInteractions.sh])
 AC_CONFIG_FILES([tests/testParser.sh], [chmod +x tests/testParser.sh])
+AC_CONFIG_FILES([tests/testSelectOutput.sh], [chmod +x tests/testSelectOutput.sh])
 
 # Save the compilation options
 AC_DEFINE_UNQUOTED([SWIFT_CONFIG_FLAGS],["$swift_config_flags"],[Flags passed to configure])
diff --git a/examples/SedovBlast_3D/sedov.yml b/examples/SedovBlast_3D/sedov.yml
index 2df2c432cef8afec49c687b643a872b06f9abb60..6b64d707a866f6abf499530105dc7af641ef7c3e 100644
--- a/examples/SedovBlast_3D/sedov.yml
+++ b/examples/SedovBlast_3D/sedov.yml
@@ -19,7 +19,7 @@ Snapshots:
   time_first:          0.    # Time of the first output (in internal units)
   delta_time:          1e-2  # Time difference between consecutive outputs (in internal units)
   compression:         4
-  
+ 
 # Parameters governing the conserved quantities statistics
 Statistics:
   delta_time:          1e-3 # Time between statistics output
@@ -33,4 +33,4 @@ SPH:
 InitialConditions:
   file_name:                    ./sedov.hdf5          
   smoothing_length_scaling:     3.33
-
+ 
\ No newline at end of file
diff --git a/examples/main.c b/examples/main.c
index 692eb53deb4d411b5cbbbd9ed8269e70aca066a1..cf30a9c2b303f843db04acb73fff2d6a70ac5f62 100644
--- a/examples/main.c
+++ b/examples/main.c
@@ -90,6 +90,8 @@ void print_help_message(void) {
   printf("  %2s %14s %s\n", "-n", "{int}",
          "Execute a fixed number of time steps. When unset use the time_end "
          "parameter to stop.");
+  printf("  %2s %14s %s\n", "-o", "{str}",
+         "Generate a default output parameter file.");
   printf("  %2s %14s %s\n", "-P", "{sec:par:val}",
          "Set parameter value and overwrites values read from the parameters "
          "file. Can be used more than once.");
@@ -195,6 +197,7 @@ int main(int argc, char *argv[]) {
   int nr_threads = 1;
   int with_verbose_timers = 0;
   int nparams = 0;
+  char output_parameters_filename[200] = "";
   char *cmdparams[PARSER_MAX_NO_OF_PARAMS];
   char paramFileName[200] = "";
   char restart_file[200] = "";
@@ -202,7 +205,7 @@ int main(int argc, char *argv[]) {
 
   /* Parse the parameters */
   int c;
-  while ((c = getopt(argc, argv, "acCdDef:FgGhMn:P:rsSt:Tv:y:Y:")) != -1)
+  while ((c = getopt(argc, argv, "acCdDef:FgGhMn:o:P:rsSt:Tv:y:Y:")) != -1)
     switch (c) {
       case 'a':
 #if defined(HAVE_SETAFFINITY) && defined(HAVE_LIBNUMA)
@@ -259,6 +262,15 @@ int main(int argc, char *argv[]) {
           return 1;
         }
         break;
+      case 'o':
+        if (sscanf(optarg, "%s", output_parameters_filename) != 1) {
+          if (myrank == 0) {
+            printf("Error parsing output fields filename");
+            print_help_message();
+          }
+          return 1;
+        }
+        break;
       case 'P':
         cmdparams[nparams] = optarg;
         nparams++;
@@ -324,6 +336,15 @@ int main(int argc, char *argv[]) {
         return 1;
         break;
     }
+
+  /* Write output parameter file */
+  if (myrank == 0 && strcmp(output_parameters_filename, "") != 0) {
+    io_write_output_field_parameter(output_parameters_filename);
+    printf("End of run.\n");
+    return 0;
+  }
+
+  /* check inputs */
   if (optind == argc - 1) {
     if (!strcpy(paramFileName, argv[optind++]))
       error("Error reading parameter file name.");
@@ -713,6 +734,9 @@ int main(int argc, char *argv[]) {
           "ICs.",
           N_total[0], N_total[2], N_total[1]);
 
+    /* Verify that the fields to dump actually exist */
+    if (myrank == 0) io_check_output_fields(params, N_total);
+
     /* Initialize the space with these data. */
     if (myrank == 0) clocks_gettime(&tic);
     space_init(&s, params, &cosmo, dim, parts, gparts, sparts, Ngas, Ngpart,
@@ -811,6 +835,7 @@ int main(int argc, char *argv[]) {
                 &cooling_func, &chemistry, &sourceterms);
     engine_config(0, &e, params, nr_nodes, myrank, nr_threads, with_aff,
                   talking, restart_file);
+
     if (myrank == 0) {
       clocks_gettime(&toc);
       message("engine_init took %.3f %s.", clocks_diff(&tic, &toc),
diff --git a/src/chemistry/EAGLE/chemistry_io.h b/src/chemistry/EAGLE/chemistry_io.h
index aab8ec240207a47289e35a711af8b245bf2b40fa..f87807579c46fe336f45e934d828318aed0377c7 100644
--- a/src/chemistry/EAGLE/chemistry_io.h
+++ b/src/chemistry/EAGLE/chemistry_io.h
@@ -30,7 +30,8 @@
  *
  * @return Returns the number of fields to read.
  */
-int chemistry_read_particles(struct part* parts, struct io_props* list) {
+INLINE static int chemistry_read_particles(struct part* parts,
+                                           struct io_props* list) {
 
   /* Nothing to read */
   return 0;
@@ -44,7 +45,8 @@ int chemistry_read_particles(struct part* parts, struct io_props* list) {
  *
  * @return Returns the number of fields to write.
  */
-int chemistry_write_particles(const struct part* parts, struct io_props* list) {
+INLINE static int chemistry_write_particles(const struct part* parts,
+                                            struct io_props* list) {
 
   /* List what we want to write */
   list[0] = io_make_output_field("ElementAbundance", FLOAT,
@@ -101,7 +103,7 @@ int chemistry_write_particles(const struct part* parts, struct io_props* list) {
  * @brief Writes the current model of SPH to the file
  * @param h_grpsph The HDF5 group in which to write
  */
-void chemistry_write_flavour(hid_t h_grp) {
+INLINE static void chemistry_write_flavour(hid_t h_grp) {
 
   io_write_attribute_s(h_grp, "Chemistry Model", "EAGLE");
   for (int elem = 0; elem < chemistry_element_count; ++elem) {
diff --git a/src/chemistry/GEAR/chemistry_io.h b/src/chemistry/GEAR/chemistry_io.h
index 0557d5c520dfc7ad5eaff2b92e6588751c072df5..2a0847bebfb8c1734f21bda2f6ad55b354a7aec9 100644
--- a/src/chemistry/GEAR/chemistry_io.h
+++ b/src/chemistry/GEAR/chemistry_io.h
@@ -48,8 +48,8 @@ chemistry_get_element_name(enum chemistry_element elem) {
  *
  * @return Returns the number of fields to read.
  */
-__attribute__((always_inline)) INLINE static int chemistry_read_particles(
-    struct part* parts, struct io_props* list) {
+INLINE static int chemistry_read_particles(struct part* parts,
+                                           struct io_props* list) {
 
   /* List what we want to read */
   list[0] = io_make_input_field(
@@ -69,13 +69,14 @@ __attribute__((always_inline)) INLINE static int chemistry_read_particles(
  *
  * @return Returns the number of fields to write.
  */
-__attribute__((always_inline)) INLINE static int chemistry_write_particles(
-    const struct part* parts, struct io_props* list) {
+INLINE static int chemistry_write_particles(const struct part* parts,
+                                            struct io_props* list) {
 
   /* List what we want to write */
   list[0] = io_make_output_field(
       "SmoothedElementAbundance", FLOAT, chemistry_element_count,
       UNIT_CONV_NO_UNITS, parts, chemistry_data.smoothed_metal_mass_fraction);
+
   list[1] = io_make_output_field("Z", FLOAT, 1, UNIT_CONV_NO_UNITS, parts,
                                  chemistry_data.Z);
 
@@ -92,8 +93,7 @@ __attribute__((always_inline)) INLINE static int chemistry_write_particles(
  * @brief Writes the current model of SPH to the file
  * @param h_grp The HDF5 group in which to write
  */
-__attribute__((always_inline)) INLINE static void chemistry_write_flavour(
-    hid_t h_grp) {
+INLINE static void chemistry_write_flavour(hid_t h_grp) {
 
   io_write_attribute_s(h_grp, "Chemistry Model", "GEAR");
   for (enum chemistry_element i = chemistry_element_O;
diff --git a/src/chemistry/none/chemistry_io.h b/src/chemistry/none/chemistry_io.h
index 142d2f75ce1487393e8689edbb9a6fdb1b1e85cd..ef7e0d8d87dfeab5978f0e86bbf6279f7901d10a 100644
--- a/src/chemistry/none/chemistry_io.h
+++ b/src/chemistry/none/chemistry_io.h
@@ -29,7 +29,8 @@
  *
  * @return Returns the number of fields to write.
  */
-int chemistry_read_particles(struct part* parts, struct io_props* list) {
+INLINE static int chemistry_read_particles(struct part* parts,
+                                           struct io_props* list) {
 
   /* update list according to hydro_io */
 
@@ -45,7 +46,8 @@ int chemistry_read_particles(struct part* parts, struct io_props* list) {
  *
  * @return Returns the number of fields to write.
  */
-int chemistry_write_particles(const struct part* parts, struct io_props* list) {
+INLINE static int chemistry_write_particles(const struct part* parts,
+                                            struct io_props* list) {
 
   /* update list according to hydro_io */
 
@@ -59,7 +61,7 @@ int chemistry_write_particles(const struct part* parts, struct io_props* list) {
  * @brief Writes the current model of SPH to the file
  * @param h_grp The HDF5 group in which to write
  */
-void chemistry_write_flavour(hid_t h_grp) {
+INLINE static void chemistry_write_flavour(hid_t h_grp) {
 
   io_write_attribute_s(h_grp, "Chemistry Model", "None");
 }
diff --git a/src/common_io.c b/src/common_io.c
index 88374e30cccfa619352aed5c7e401bf163d6db63..d5aaa71ca2b9d8ed157559d5edc1f406df454e86 100644
--- a/src/common_io.c
+++ b/src/common_io.c
@@ -25,12 +25,17 @@
 #include "common_io.h"
 
 /* Local includes. */
+#include "chemistry_io.h"
 #include "engine.h"
 #include "error.h"
+#include "gravity_io.h"
 #include "hydro.h"
+#include "hydro_io.h"
 #include "io_properties.h"
 #include "kernel_hydro.h"
 #include "part.h"
+#include "part_type.h"
+#include "stars_io.h"
 #include "threadpool.h"
 #include "units.h"
 #include "version.h"
@@ -806,3 +811,157 @@ void io_collect_dm_gparts(const struct gpart* const gparts, size_t Ntot,
     error("Collected the wrong number of dm particles (%zu vs. %zu expected)",
           count, Ndm);
 }
+
+/**
+ * @brief Verify the io parameter file
+ *
+ * @param params The #swift_params
+ * @param N_total The total number of each particle type.
+ */
+void io_check_output_fields(const struct swift_params* params,
+                            const long long N_total[3]) {
+
+  /* Create some fake particles as arguments for the writing routines */
+  struct part p;
+  struct xpart xp;
+  struct spart sp;
+  struct gpart gp;
+
+  /* Copy N_total to array with length == 6 */
+  const long long nr_total[swift_type_count] = {
+      N_total[0], N_total[1], N_total[2], 0, 0, 0};
+
+  /* Loop over all particle types to check the fields */
+  for (int ptype = 0; ptype < swift_type_count; ptype++) {
+
+    int num_fields = 0;
+    struct io_props list[100];
+
+    /* Don't do anything if no particle of this kind */
+    if (nr_total[ptype] == 0) continue;
+
+    /* Gather particle fields from the particle structures */
+    switch (ptype) {
+
+      case swift_type_gas:
+        hydro_write_particles(&p, &xp, list, &num_fields);
+        num_fields += chemistry_write_particles(&p, list + num_fields);
+        break;
+
+      case swift_type_dark_matter:
+        darkmatter_write_particles(&gp, list, &num_fields);
+        break;
+
+      case swift_type_star:
+        star_write_particles(&sp, list, &num_fields);
+        break;
+
+      default:
+        error("Particle Type %d not yet supported. Aborting", ptype);
+    }
+
+    /* loop over each parameter */
+    for (int param_id = 0; param_id < params->paramCount; param_id++) {
+      const char* param_name = params->data[param_id].name;
+
+      char section_name[PARSER_MAX_LINE_SIZE];
+
+      /* Skip if wrong section */
+      sprintf(section_name, "SelectOutput:");
+      if (strstr(param_name, section_name) == NULL) continue;
+
+      /* Skip if wrong particle type */
+      sprintf(section_name, "_%s", part_type_names[ptype]);
+      if (strstr(param_name, section_name) == NULL) continue;
+
+      int found = 0;
+
+      /* loop over each possible output field */
+      for (int field_id = 0; field_id < num_fields; field_id++) {
+        char field_name[PARSER_MAX_LINE_SIZE];
+        sprintf(field_name, "SelectOutput:%s_%s", list[field_id].name,
+                part_type_names[ptype]);
+
+        if (strcmp(param_name, field_name) == 0) {
+          found = 1;
+          /* check if correct input */
+          int retParam = 0;
+          char str[PARSER_MAX_LINE_SIZE];
+          sscanf(params->data[param_id].value, "%d%s", &retParam, str);
+
+          /* Check that we have a 0 or 1 */
+          if (retParam != 0 && retParam != 1)
+            message(
+                "WARNING: Unexpected input for %s. Received %i but expect 0 or "
+                "1. ",
+                field_name, retParam);
+
+          /* Found it, so move to the next one. */
+          break;
+        }
+      }
+      if (!found)
+        message(
+            "WARNING: Trying to dump particle field '%s' (read from '%s') that "
+            "does not exist.",
+            param_name, params->fileName);
+    }
+  }
+}
+
+/**
+ * @brief Write the output field parameters file
+ *
+ * @param filename The file to write
+ */
+void io_write_output_field_parameter(const char* filename) {
+
+  FILE* file = fopen(filename, "w");
+  if (file == NULL) error("Error opening file '%s'", filename);
+
+  /* Loop over all particle types */
+  fprintf(file, "SelectOutput:\n");
+  for (int ptype = 0; ptype < swift_type_count; ptype++) {
+
+    int num_fields = 0;
+    struct io_props list[100];
+
+    /* Write particle fields from the particle structure */
+    switch (ptype) {
+
+      case swift_type_gas:
+        hydro_write_particles(NULL, NULL, list, &num_fields);
+        num_fields += chemistry_write_particles(NULL, list + num_fields);
+        break;
+
+      case swift_type_dark_matter:
+        darkmatter_write_particles(NULL, list, &num_fields);
+        break;
+
+      case swift_type_star:
+        star_write_particles(NULL, list, &num_fields);
+        break;
+
+      default:
+        break;
+    }
+
+    if (num_fields == 0) continue;
+
+    /* Output a header for that particle type */
+    fprintf(file, "  # Particle Type %s\n", part_type_names[ptype]);
+
+    /* Write all the fields of this particle type */
+    for (int i = 0; i < num_fields; ++i)
+      fprintf(file, "  %s_%s: 1\n", list[i].name, part_type_names[ptype]);
+
+    fprintf(file, "\n");
+  }
+
+  fclose(file);
+
+  printf(
+      "List of valid ouput fields for the particle in snapshots dumped in "
+      "'%s'.\n",
+      filename);
+}
diff --git a/src/common_io.h b/src/common_io.h
index d9e676db934b58ee476f18894acf55c4d38344f9..f26a635a66f40424984238e586fcdf5bc752fc99 100644
--- a/src/common_io.h
+++ b/src/common_io.h
@@ -100,4 +100,9 @@ void io_duplicate_star_gparts(struct threadpool* tp, struct spart* const sparts,
                               struct gpart* const gparts, size_t Nstars,
                               size_t Ndm);
 
+void io_check_output_fields(const struct swift_params* params,
+                            const long long N_total[3]);
+
+void io_write_output_field_parameter(const char* filename);
+
 #endif /* SWIFT_COMMON_IO_H */
diff --git a/src/gravity/Default/gravity_io.h b/src/gravity/Default/gravity_io.h
index 7b8ec2c8fef04cc4a4fc6836d8bb895b24d3c41f..5016803b1fb1e4c6482cde91fae555d4d3e28c9d 100644
--- a/src/gravity/Default/gravity_io.h
+++ b/src/gravity/Default/gravity_io.h
@@ -21,8 +21,8 @@
 
 #include "io_properties.h"
 
-void convert_gpart_pos(const struct engine* e, const struct gpart* gp,
-                       double* ret) {
+INLINE static void convert_gpart_pos(const struct engine* e,
+                                     const struct gpart* gp, double* ret) {
 
   if (e->s->periodic) {
     ret[0] = box_wrap(gp->x[0], 0.0, e->s->dim[0]);
@@ -35,8 +35,8 @@ void convert_gpart_pos(const struct engine* e, const struct gpart* gp,
   }
 }
 
-void convert_gpart_vel(const struct engine* e, const struct gpart* gp,
-                       float* ret) {
+INLINE static void convert_gpart_vel(const struct engine* e,
+                                     const struct gpart* gp, float* ret) {
 
   const int with_cosmology = (e->policy & engine_policy_cosmology);
   const struct cosmology* cosmo = e->cosmology;
@@ -74,8 +74,9 @@ void convert_gpart_vel(const struct engine* e, const struct gpart* gp,
  * @param list The list of i/o properties to read.
  * @param num_fields The number of i/o fields to read.
  */
-void darkmatter_read_particles(struct gpart* gparts, struct io_props* list,
-                               int* num_fields) {
+INLINE static void darkmatter_read_particles(struct gpart* gparts,
+                                             struct io_props* list,
+                                             int* num_fields) {
 
   /* Say how much we want to read */
   *num_fields = 4;
@@ -98,8 +99,9 @@ void darkmatter_read_particles(struct gpart* gparts, struct io_props* list,
  * @param list The list of i/o properties to write.
  * @param num_fields The number of i/o fields to write.
  */
-void darkmatter_write_particles(const struct gpart* gparts,
-                                struct io_props* list, int* num_fields) {
+INLINE static void darkmatter_write_particles(const struct gpart* gparts,
+                                              struct io_props* list,
+                                              int* num_fields) {
 
   /* Say how much we want to write */
   *num_fields = 5;
diff --git a/src/hydro/Default/hydro_io.h b/src/hydro/Default/hydro_io.h
index c7a0953eba5e5c713caa4f01cb9e8bd20f6b065f..d6b6f36ba23cee02690091611fb485ad54b47b4d 100644
--- a/src/hydro/Default/hydro_io.h
+++ b/src/hydro/Default/hydro_io.h
@@ -31,8 +31,9 @@
  * @param list The list of i/o properties to read.
  * @param num_fields The number of i/o fields to read.
  */
-void hydro_read_particles(struct part* parts, struct io_props* list,
-                          int* num_fields) {
+INLINE static void hydro_read_particles(struct part* parts,
+                                        struct io_props* list,
+                                        int* num_fields) {
 
   *num_fields = 8;
 
@@ -55,8 +56,9 @@ void hydro_read_particles(struct part* parts, struct io_props* list,
                                 UNIT_CONV_DENSITY, parts, rho);
 }
 
-void convert_part_pos(const struct engine* e, const struct part* p,
-                      const struct xpart* xp, double* ret) {
+INLINE static void convert_part_pos(const struct engine* e,
+                                    const struct part* p,
+                                    const struct xpart* xp, double* ret) {
 
   if (e->s->periodic) {
     ret[0] = box_wrap(p->x[0], 0.0, e->s->dim[0]);
@@ -69,8 +71,9 @@ void convert_part_pos(const struct engine* e, const struct part* p,
   }
 }
 
-void convert_part_vel(const struct engine* e, const struct part* p,
-                      const struct xpart* xp, float* ret) {
+INLINE static void convert_part_vel(const struct engine* e,
+                                    const struct part* p,
+                                    const struct xpart* xp, float* ret) {
 
   const int with_cosmology = (e->policy & engine_policy_cosmology);
   const struct cosmology* cosmo = e->cosmology;
@@ -103,8 +106,9 @@ void convert_part_vel(const struct engine* e, const struct part* p,
   ret[2] *= cosmo->a2_inv;
 }
 
-void convert_part_potential(const struct engine* e, const struct part* p,
-                            const struct xpart* xp, float* ret) {
+INLINE static void convert_part_potential(const struct engine* e,
+                                          const struct part* p,
+                                          const struct xpart* xp, float* ret) {
 
   if (p->gpart != NULL)
     ret[0] = gravity_get_comoving_potential(p->gpart);
@@ -119,8 +123,10 @@ void convert_part_potential(const struct engine* e, const struct part* p,
  * @param list The list of i/o properties to write.
  * @param num_fields The number of i/o fields to write.
  */
-void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
-                           struct io_props* list, int* num_fields) {
+INLINE static void hydro_write_particles(const struct part* parts,
+                                         const struct xpart* xparts,
+                                         struct io_props* list,
+                                         int* num_fields) {
 
   *num_fields = 8;
 
@@ -149,7 +155,7 @@ void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
  * @brief Writes the current model of SPH to the file
  * @param h_grpsph The HDF5 group in which to write
  */
-void hydro_write_flavour(hid_t h_grpsph) {
+INLINE static void hydro_write_flavour(hid_t h_grpsph) {
 
   /* Viscosity and thermal conduction */
   io_write_attribute_s(h_grpsph, "Thermal Conductivity Model",
@@ -178,6 +184,6 @@ void hydro_write_flavour(hid_t h_grpsph) {
  *
  * @return 1 if entropy is in 'internal energy', 0 otherwise.
  */
-int writeEntropyFlag(void) { return 0; }
+INLINE static int writeEntropyFlag(void) { return 0; }
 
 #endif /* SWIFT_DEFAULT_HYDRO_IO_H */
diff --git a/src/hydro/Gadget2/hydro_io.h b/src/hydro/Gadget2/hydro_io.h
index b6e2c32230cac682c7ca907dcada971117be55c0..8f503230e1b2014c506e4c476f8597b5b77f85a1 100644
--- a/src/hydro/Gadget2/hydro_io.h
+++ b/src/hydro/Gadget2/hydro_io.h
@@ -31,8 +31,9 @@
  * @param list The list of i/o properties to read.
  * @param num_fields The number of i/o fields to read.
  */
-void hydro_read_particles(struct part* parts, struct io_props* list,
-                          int* num_fields) {
+INLINE static void hydro_read_particles(struct part* parts,
+                                        struct io_props* list,
+                                        int* num_fields) {
 
   *num_fields = 8;
 
@@ -55,20 +56,21 @@ void hydro_read_particles(struct part* parts, struct io_props* list,
                                 UNIT_CONV_DENSITY, parts, rho);
 }
 
-void convert_part_u(const struct engine* e, const struct part* p,
-                    const struct xpart* xp, float* ret) {
+INLINE static void convert_part_u(const struct engine* e, const struct part* p,
+                                  const struct xpart* xp, float* ret) {
 
   ret[0] = hydro_get_comoving_internal_energy(p);
 }
 
-void convert_part_P(const struct engine* e, const struct part* p,
-                    const struct xpart* xp, float* ret) {
+INLINE static void convert_part_P(const struct engine* e, const struct part* p,
+                                  const struct xpart* xp, float* ret) {
 
   ret[0] = hydro_get_comoving_pressure(p);
 }
 
-void convert_part_pos(const struct engine* e, const struct part* p,
-                      const struct xpart* xp, double* ret) {
+INLINE static void convert_part_pos(const struct engine* e,
+                                    const struct part* p,
+                                    const struct xpart* xp, double* ret) {
 
   if (e->s->periodic) {
     ret[0] = box_wrap(p->x[0], 0.0, e->s->dim[0]);
@@ -81,8 +83,9 @@ void convert_part_pos(const struct engine* e, const struct part* p,
   }
 }
 
-void convert_part_vel(const struct engine* e, const struct part* p,
-                      const struct xpart* xp, float* ret) {
+INLINE static void convert_part_vel(const struct engine* e,
+                                    const struct part* p,
+                                    const struct xpart* xp, float* ret) {
 
   const int with_cosmology = (e->policy & engine_policy_cosmology);
   const struct cosmology* cosmo = e->cosmology;
@@ -115,8 +118,9 @@ void convert_part_vel(const struct engine* e, const struct part* p,
   ret[2] *= cosmo->a2_inv;
 }
 
-void convert_part_potential(const struct engine* e, const struct part* p,
-                            const struct xpart* xp, float* ret) {
+INLINE static void convert_part_potential(const struct engine* e,
+                                          const struct part* p,
+                                          const struct xpart* xp, float* ret) {
 
   if (p->gpart != NULL)
     ret[0] = gravity_get_comoving_potential(p->gpart);
@@ -131,8 +135,10 @@ void convert_part_potential(const struct engine* e, const struct part* p,
  * @param list The list of i/o properties to write.
  * @param num_fields The number of i/o fields to write.
  */
-void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
-                           struct io_props* list, int* num_fields) {
+INLINE static void hydro_write_particles(const struct part* parts,
+                                         const struct xpart* xparts,
+                                         struct io_props* list,
+                                         int* num_fields) {
 
   *num_fields = 10;
 
@@ -185,7 +191,7 @@ void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
  * @brief Writes the current model of SPH to the file
  * @param h_grpsph The HDF5 group in which to write
  */
-void hydro_write_flavour(hid_t h_grpsph) {
+INLINE static void hydro_write_flavour(hid_t h_grpsph) {
 
   /* Viscosity and thermal conduction */
   io_write_attribute_s(h_grpsph, "Thermal Conductivity Model",
@@ -202,6 +208,6 @@ void hydro_write_flavour(hid_t h_grpsph) {
  *
  * @return 1 if entropy is in 'internal energy', 0 otherwise.
  */
-int writeEntropyFlag(void) { return 0; }
+INLINE static int writeEntropyFlag(void) { return 0; }
 
 #endif /* SWIFT_GADGET2_HYDRO_IO_H */
diff --git a/src/hydro/GizmoMFM/hydro_io.h b/src/hydro/GizmoMFM/hydro_io.h
index 94127b3fa954446e8d8de804430ca7319ff35381..6ce5a23eeed66a1f6a7205d4683ecaeb91e05259 100644
--- a/src/hydro/GizmoMFM/hydro_io.h
+++ b/src/hydro/GizmoMFM/hydro_io.h
@@ -40,8 +40,9 @@
  * @param list The list of i/o properties to read.
  * @param num_fields The number of i/o fields to read.
  */
-void hydro_read_particles(struct part* parts, struct io_props* list,
-                          int* num_fields) {
+INLINE static void hydro_read_particles(struct part* parts,
+                                        struct io_props* list,
+                                        int* num_fields) {
 
   *num_fields = 8;
 
@@ -72,8 +73,8 @@ void hydro_read_particles(struct part* parts, struct io_props* list,
  * @param p Particle.
  * @param ret (return) Internal energy of the particle
  */
-void convert_u(const struct engine* e, const struct part* p,
-               const struct xpart* xp, float* ret) {
+INLINE static void convert_u(const struct engine* e, const struct part* p,
+                             const struct xpart* xp, float* ret) {
 
   ret[0] = hydro_get_comoving_internal_energy(p);
 }
@@ -85,8 +86,8 @@ void convert_u(const struct engine* e, const struct part* p,
  * @param p Particle.
  * @param ret (return) Entropic function of the particle
  */
-void convert_A(const struct engine* e, const struct part* p,
-               const struct xpart* xp, float* ret) {
+INLINE static void convert_A(const struct engine* e, const struct part* p,
+                             const struct xpart* xp, float* ret) {
   ret[0] = hydro_get_comoving_entropy(p);
 }
 
@@ -97,8 +98,8 @@ void convert_A(const struct engine* e, const struct part* p,
  * @param p Particle.
  * @return Total energy of the particle
  */
-void convert_Etot(const struct engine* e, const struct part* p,
-                  const struct xpart* xp, float* ret) {
+INLINE static void convert_Etot(const struct engine* e, const struct part* p,
+                                const struct xpart* xp, float* ret) {
 #ifdef GIZMO_TOTAL_ENERGY
   ret[0] = p->conserved.energy;
 #else
@@ -112,8 +113,9 @@ void convert_Etot(const struct engine* e, const struct part* p,
 #endif
 }
 
-void convert_part_pos(const struct engine* e, const struct part* p,
-                      const struct xpart* xp, double* ret) {
+INLINE static void convert_part_pos(const struct engine* e,
+                                    const struct part* p,
+                                    const struct xpart* xp, double* ret) {
 
   if (e->s->periodic) {
     ret[0] = box_wrap(p->x[0], 0.0, e->s->dim[0]);
@@ -126,8 +128,9 @@ void convert_part_pos(const struct engine* e, const struct part* p,
   }
 }
 
-void convert_part_vel(const struct engine* e, const struct part* p,
-                      const struct xpart* xp, float* ret) {
+INLINE static void convert_part_vel(const struct engine* e,
+                                    const struct part* p,
+                                    const struct xpart* xp, float* ret) {
 
   const int with_cosmology = (e->policy & engine_policy_cosmology);
   const struct cosmology* cosmo = e->cosmology;
@@ -160,8 +163,9 @@ void convert_part_vel(const struct engine* e, const struct part* p,
   ret[2] *= cosmo->a2_inv;
 }
 
-void convert_part_potential(const struct engine* e, const struct part* p,
-                            const struct xpart* xp, float* ret) {
+INLINE static void convert_part_potential(const struct engine* e,
+                                          const struct part* p,
+                                          const struct xpart* xp, float* ret) {
 
   if (p->gpart != NULL)
     ret[0] = gravity_get_comoving_potential(p->gpart);
@@ -176,8 +180,10 @@ void convert_part_potential(const struct engine* e, const struct part* p,
  * @param list The list of i/o properties to write.
  * @param num_fields The number of i/o fields to write.
  */
-void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
-                           struct io_props* list, int* num_fields) {
+INLINE static void hydro_write_particles(const struct part* parts,
+                                         const struct xpart* xparts,
+                                         struct io_props* list,
+                                         int* num_fields) {
 
   *num_fields = 11;
 
@@ -215,7 +221,7 @@ void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
  * @brief Writes the current model of SPH to the file
  * @param h_grpsph The HDF5 group in which to write
  */
-void hydro_write_flavour(hid_t h_grpsph) {
+INLINE static void hydro_write_flavour(hid_t h_grpsph) {
   /* Gradient information */
   io_write_attribute_s(h_grpsph, "Gradient reconstruction model",
                        HYDRO_GRADIENT_IMPLEMENTATION);
@@ -239,6 +245,6 @@ void hydro_write_flavour(hid_t h_grpsph) {
  *
  * @return 1 if entropy is in 'internal energy', 0 otherwise.
  */
-int writeEntropyFlag(void) { return 0; }
+INLINE static int writeEntropyFlag(void) { return 0; }
 
 #endif /* SWIFT_GIZMO_MFM_HYDRO_IO_H */
diff --git a/src/hydro/GizmoMFV/hydro_io.h b/src/hydro/GizmoMFV/hydro_io.h
index 18f023acbaad7bddbc4b1531d47987b5faf3bac9..7b1c651698bafad2740e179f77e4f5c3c7d1f677 100644
--- a/src/hydro/GizmoMFV/hydro_io.h
+++ b/src/hydro/GizmoMFV/hydro_io.h
@@ -40,8 +40,9 @@
  * @param list The list of i/o properties to read.
  * @param num_fields The number of i/o fields to read.
  */
-void hydro_read_particles(struct part* parts, struct io_props* list,
-                          int* num_fields) {
+INLINE static void hydro_read_particles(struct part* parts,
+                                        struct io_props* list,
+                                        int* num_fields) {
 
   *num_fields = 8;
 
@@ -72,8 +73,8 @@ void hydro_read_particles(struct part* parts, struct io_props* list,
  * @param p Particle.
  * @param ret (return) Internal energy of the particle
  */
-void convert_u(const struct engine* e, const struct part* p,
-               const struct xpart* xp, float* ret) {
+INLINE static void convert_u(const struct engine* e, const struct part* p,
+                             const struct xpart* xp, float* ret) {
 
   ret[0] = hydro_get_comoving_internal_energy(p);
 }
@@ -85,8 +86,8 @@ void convert_u(const struct engine* e, const struct part* p,
  * @param p Particle.
  * @param ret (return) Entropic function of the particle
  */
-void convert_A(const struct engine* e, const struct part* p,
-               const struct xpart* xp, float* ret) {
+INLINE static void convert_A(const struct engine* e, const struct part* p,
+                             const struct xpart* xp, float* ret) {
   ret[0] = hydro_get_comoving_entropy(p);
 }
 
@@ -97,8 +98,8 @@ void convert_A(const struct engine* e, const struct part* p,
  * @param p Particle.
  * @return Total energy of the particle
  */
-void convert_Etot(const struct engine* e, const struct part* p,
-                  const struct xpart* xp, float* ret) {
+INLINE static void convert_Etot(const struct engine* e, const struct part* p,
+                                const struct xpart* xp, float* ret) {
 #ifdef GIZMO_TOTAL_ENERGY
   ret[0] = p->conserved.energy;
 #else
@@ -112,8 +113,9 @@ void convert_Etot(const struct engine* e, const struct part* p,
 #endif
 }
 
-void convert_part_pos(const struct engine* e, const struct part* p,
-                      const struct xpart* xp, double* ret) {
+INLINE static void convert_part_pos(const struct engine* e,
+                                    const struct part* p,
+                                    const struct xpart* xp, double* ret) {
 
   if (e->s->periodic) {
     ret[0] = box_wrap(p->x[0], 0.0, e->s->dim[0]);
@@ -126,8 +128,9 @@ void convert_part_pos(const struct engine* e, const struct part* p,
   }
 }
 
-void convert_part_vel(const struct engine* e, const struct part* p,
-                      const struct xpart* xp, float* ret) {
+INLINE static void convert_part_vel(const struct engine* e,
+                                    const struct part* p,
+                                    const struct xpart* xp, float* ret) {
 
   const int with_cosmology = (e->policy & engine_policy_cosmology);
   const struct cosmology* cosmo = e->cosmology;
@@ -160,8 +163,9 @@ void convert_part_vel(const struct engine* e, const struct part* p,
   ret[2] *= cosmo->a2_inv;
 }
 
-void convert_part_potential(const struct engine* e, const struct part* p,
-                            const struct xpart* xp, float* ret) {
+INLINE static void convert_part_potential(const struct engine* e,
+                                          const struct part* p,
+                                          const struct xpart* xp, float* ret) {
 
   if (p->gpart != NULL)
     ret[0] = gravity_get_comoving_potential(p->gpart);
@@ -176,8 +180,10 @@ void convert_part_potential(const struct engine* e, const struct part* p,
  * @param list The list of i/o properties to write.
  * @param num_fields The number of i/o fields to write.
  */
-void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
-                           struct io_props* list, int* num_fields) {
+INLINE static void hydro_write_particles(const struct part* parts,
+                                         const struct xpart* xparts,
+                                         struct io_props* list,
+                                         int* num_fields) {
 
   *num_fields = 11;
 
@@ -215,7 +221,7 @@ void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
  * @brief Writes the current model of SPH to the file
  * @param h_grpsph The HDF5 group in which to write
  */
-void hydro_write_flavour(hid_t h_grpsph) {
+INLINE static void hydro_write_flavour(hid_t h_grpsph) {
   /* Gradient information */
   io_write_attribute_s(h_grpsph, "Gradient reconstruction model",
                        HYDRO_GRADIENT_IMPLEMENTATION);
@@ -239,6 +245,6 @@ void hydro_write_flavour(hid_t h_grpsph) {
  *
  * @return 1 if entropy is in 'internal energy', 0 otherwise.
  */
-int writeEntropyFlag(void) { return 0; }
+INLINE static int writeEntropyFlag(void) { return 0; }
 
 #endif /* SWIFT_GIZMO_MFV_HYDRO_IO_H */
diff --git a/src/hydro/Minimal/hydro_io.h b/src/hydro/Minimal/hydro_io.h
index 846462a287ca9f90b762c55a5c0d12d88632d7ff..833103e7ab00738070cfe3ef8fc4ce4a2b8d6348 100644
--- a/src/hydro/Minimal/hydro_io.h
+++ b/src/hydro/Minimal/hydro_io.h
@@ -45,8 +45,9 @@
  * @param list The list of i/o properties to read.
  * @param num_fields The number of i/o fields to read.
  */
-void hydro_read_particles(struct part* parts, struct io_props* list,
-                          int* num_fields) {
+INLINE static void hydro_read_particles(struct part* parts,
+                                        struct io_props* list,
+                                        int* num_fields) {
 
   *num_fields = 8;
 
@@ -69,20 +70,21 @@ void hydro_read_particles(struct part* parts, struct io_props* list,
                                 UNIT_CONV_DENSITY, parts, rho);
 }
 
-void convert_S(const struct engine* e, const struct part* p,
-               const struct xpart* xp, float* ret) {
+INLINE static void convert_S(const struct engine* e, const struct part* p,
+                             const struct xpart* xp, float* ret) {
 
   ret[0] = hydro_get_comoving_entropy(p);
 }
 
-void convert_P(const struct engine* e, const struct part* p,
-               const struct xpart* xp, float* ret) {
+INLINE static void convert_P(const struct engine* e, const struct part* p,
+                             const struct xpart* xp, float* ret) {
 
   ret[0] = hydro_get_comoving_pressure(p);
 }
 
-void convert_part_pos(const struct engine* e, const struct part* p,
-                      const struct xpart* xp, double* ret) {
+INLINE static void convert_part_pos(const struct engine* e,
+                                    const struct part* p,
+                                    const struct xpart* xp, double* ret) {
 
   if (e->s->periodic) {
     ret[0] = box_wrap(p->x[0], 0.0, e->s->dim[0]);
@@ -95,8 +97,9 @@ void convert_part_pos(const struct engine* e, const struct part* p,
   }
 }
 
-void convert_part_vel(const struct engine* e, const struct part* p,
-                      const struct xpart* xp, float* ret) {
+INLINE static void convert_part_vel(const struct engine* e,
+                                    const struct part* p,
+                                    const struct xpart* xp, float* ret) {
 
   const int with_cosmology = (e->policy & engine_policy_cosmology);
   const struct cosmology* cosmo = e->cosmology;
@@ -129,8 +132,9 @@ void convert_part_vel(const struct engine* e, const struct part* p,
   ret[2] *= cosmo->a2_inv;
 }
 
-void convert_part_potential(const struct engine* e, const struct part* p,
-                            const struct xpart* xp, float* ret) {
+INLINE static void convert_part_potential(const struct engine* e,
+                                          const struct part* p,
+                                          const struct xpart* xp, float* ret) {
 
   if (p->gpart != NULL)
     ret[0] = gravity_get_comoving_potential(p->gpart);
@@ -146,8 +150,10 @@ void convert_part_potential(const struct engine* e, const struct part* p,
  * @param list The list of i/o properties to write.
  * @param num_fields The number of i/o fields to write.
  */
-void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
-                           struct io_props* list, int* num_fields) {
+INLINE static void hydro_write_particles(const struct part* parts,
+                                         const struct xpart* xparts,
+                                         struct io_props* list,
+                                         int* num_fields) {
 
   *num_fields = 10;
 
@@ -182,7 +188,7 @@ void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
  * @brief Writes the current model of SPH to the file
  * @param h_grpsph The HDF5 group in which to write
  */
-void hydro_write_flavour(hid_t h_grpsph) {
+INLINE static void hydro_write_flavour(hid_t h_grpsph) {
 
   /* Viscosity and thermal conduction */
   /* Nothing in this minimal model... */
@@ -200,6 +206,6 @@ void hydro_write_flavour(hid_t h_grpsph) {
  *
  * @return 1 if entropy is in 'internal energy', 0 otherwise.
  */
-int writeEntropyFlag(void) { return 0; }
+INLINE static int writeEntropyFlag(void) { return 0; }
 
 #endif /* SWIFT_MINIMAL_HYDRO_IO_H */
diff --git a/src/hydro/MinimalMultiMat/hydro_io.h b/src/hydro/MinimalMultiMat/hydro_io.h
index a1b6afa785aee3d5e0bfd8043e259ff65ccfcc70..5542b93ae43d88c93ed2e2ebddc1a6baef82a557 100644
--- a/src/hydro/MinimalMultiMat/hydro_io.h
+++ b/src/hydro/MinimalMultiMat/hydro_io.h
@@ -46,8 +46,9 @@
  * @param list The list of i/o properties to read.
  * @param num_fields The number of i/o fields to read.
  */
-void hydro_read_particles(struct part* parts, struct io_props* list,
-                          int* num_fields) {
+INLINE static void hydro_read_particles(struct part* parts,
+                                        struct io_props* list,
+                                        int* num_fields) {
 
   *num_fields = 9;
 
@@ -72,20 +73,21 @@ void hydro_read_particles(struct part* parts, struct io_props* list,
       io_make_input_field("MaterialID", INT, 1, OPTIONAL, 1, parts, mat_id);
 }
 
-void convert_S(const struct engine* e, const struct part* p,
-               const struct xpart* xp, float* ret) {
+INLINE static void convert_S(const struct engine* e, const struct part* p,
+                             const struct xpart* xp, float* ret) {
 
   ret[0] = hydro_get_comoving_entropy(p);
 }
 
-void convert_P(const struct engine* e, const struct part* p,
-               const struct xpart* xp, float* ret) {
+INLINE static void convert_P(const struct engine* e, const struct part* p,
+                             const struct xpart* xp, float* ret) {
 
   ret[0] = hydro_get_comoving_pressure(p);
 }
 
-void convert_part_pos(const struct engine* e, const struct part* p,
-                      const struct xpart* xp, double* ret) {
+INLINE static void convert_part_pos(const struct engine* e,
+                                    const struct part* p,
+                                    const struct xpart* xp, double* ret) {
 
   if (e->s->periodic) {
     ret[0] = box_wrap(p->x[0], 0.0, e->s->dim[0]);
@@ -98,8 +100,9 @@ void convert_part_pos(const struct engine* e, const struct part* p,
   }
 }
 
-void convert_part_vel(const struct engine* e, const struct part* p,
-                      const struct xpart* xp, float* ret) {
+INLINE static void convert_part_vel(const struct engine* e,
+                                    const struct part* p,
+                                    const struct xpart* xp, float* ret) {
 
   const int with_cosmology = (e->policy & engine_policy_cosmology);
   const struct cosmology* cosmo = e->cosmology;
@@ -132,8 +135,9 @@ void convert_part_vel(const struct engine* e, const struct part* p,
   ret[2] *= cosmo->a2_inv;
 }
 
-void convert_part_potential(const struct engine* e, const struct part* p,
-                            const struct xpart* xp, float* ret) {
+INLINE static void convert_part_potential(const struct engine* e,
+                                          const struct part* p,
+                                          const struct xpart* xp, float* ret) {
 
   if (p->gpart != NULL)
     ret[0] = gravity_get_comoving_potential(p->gpart);
@@ -149,8 +153,10 @@ void convert_part_potential(const struct engine* e, const struct part* p,
  * @param list The list of i/o properties to write.
  * @param num_fields The number of i/o fields to write.
  */
-void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
-                           struct io_props* list, int* num_fields) {
+INLINE static void hydro_write_particles(const struct part* parts,
+                                         const struct xpart* xparts,
+                                         struct io_props* list,
+                                         int* num_fields) {
 
   *num_fields = 11;
 
@@ -186,7 +192,7 @@ void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
  * @brief Writes the current model of SPH to the file
  * @param h_grpsph The HDF5 group in which to write
  */
-void hydro_write_flavour(hid_t h_grpsph) {
+INLINE static void hydro_write_flavour(hid_t h_grpsph) {
 
   /* Viscosity and thermal conduction */
   /* Nothing in this minimal model... */
@@ -204,6 +210,6 @@ void hydro_write_flavour(hid_t h_grpsph) {
  *
  * @return 1 if entropy is in 'internal energy', 0 otherwise.
  */
-int writeEntropyFlag(void) { return 0; }
+INLINE static int writeEntropyFlag(void) { return 0; }
 
 #endif /* SWIFT_MINIMAL_MULTI_MAT_HYDRO_IO_H */
diff --git a/src/hydro/PressureEnergy/hydro_io.h b/src/hydro/PressureEnergy/hydro_io.h
index 521f0ba28efc18716e6fb46c24a9d86bbf6d9d5a..4f168ed81961438b996f7973c9b04a8525ecf039 100644
--- a/src/hydro/PressureEnergy/hydro_io.h
+++ b/src/hydro/PressureEnergy/hydro_io.h
@@ -43,8 +43,9 @@
  * @param list The list of i/o properties to read.
  * @param num_fields The number of i/o fields to read.
  */
-void hydro_read_particles(struct part* parts, struct io_props* list,
-                          int* num_fields) {
+INLINE static void hydro_read_particles(struct part* parts,
+                                        struct io_props* list,
+                                        int* num_fields) {
 
   *num_fields = 8;
 
@@ -67,26 +68,27 @@ void hydro_read_particles(struct part* parts, struct io_props* list,
                                 UNIT_CONV_DENSITY, parts, rho);
 }
 
-void convert_u(const struct engine* e, const struct part* p,
-               const struct xpart* xp, float* ret) {
+INLINE static void convert_u(const struct engine* e, const struct part* p,
+                             const struct xpart* xp, float* ret) {
 
   ret[0] = hydro_get_comoving_internal_energy(p);
 }
 
-void convert_S(const struct engine* e, const struct part* p,
-               const struct xpart* xp, float* ret) {
+INLINE static void convert_S(const struct engine* e, const struct part* p,
+                             const struct xpart* xp, float* ret) {
 
   ret[0] = hydro_get_comoving_entropy(p);
 }
 
-void convert_P(const struct engine* e, const struct part* p,
-               const struct xpart* xp, float* ret) {
+INLINE static void convert_P(const struct engine* e, const struct part* p,
+                             const struct xpart* xp, float* ret) {
 
   ret[0] = hydro_get_comoving_pressure(p);
 }
 
-void convert_part_pos(const struct engine* e, const struct part* p,
-                      const struct xpart* xp, double* ret) {
+INLINE static void convert_part_pos(const struct engine* e,
+                                    const struct part* p,
+                                    const struct xpart* xp, double* ret) {
 
   if (e->s->periodic) {
     ret[0] = box_wrap(p->x[0], 0.0, e->s->dim[0]);
@@ -99,8 +101,9 @@ void convert_part_pos(const struct engine* e, const struct part* p,
   }
 }
 
-void convert_part_vel(const struct engine* e, const struct part* p,
-                      const struct xpart* xp, float* ret) {
+INLINE static void convert_part_vel(const struct engine* e,
+                                    const struct part* p,
+                                    const struct xpart* xp, float* ret) {
 
   const int with_cosmology = (e->policy & engine_policy_cosmology);
   const struct cosmology* cosmo = e->cosmology;
@@ -140,8 +143,10 @@ void convert_part_vel(const struct engine* e, const struct part* p,
  * @param list The list of i/o properties to write.
  * @param num_fields The number of i/o fields to write.
  */
-void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
-                           struct io_props* list, int* num_fields) {
+INLINE static void hydro_write_particles(const struct part* parts,
+                                         const struct xpart* xparts,
+                                         struct io_props* list,
+                                         int* num_fields) {
 
   *num_fields = 9;
 
@@ -173,7 +178,7 @@ void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
  * @brief Writes the current model of SPH to the file
  * @param h_grpsph The HDF5 group in which to write
  */
-void hydro_write_flavour(hid_t h_grpsph) {
+INLINE static void hydro_write_flavour(hid_t h_grpsph) {
 
   /* Viscosity and thermal conduction */
   /* Nothing in this minimal model... */
@@ -191,6 +196,6 @@ void hydro_write_flavour(hid_t h_grpsph) {
  *
  * @return 1 if entropy is in 'internal energy', 0 otherwise.
  */
-int writeEntropyFlag(void) { return 0; }
+INLINE static int writeEntropyFlag(void) { return 0; }
 
 #endif /* SWIFT_MINIMAL_HYDRO_IO_H */
diff --git a/src/hydro/PressureEntropy/hydro_io.h b/src/hydro/PressureEntropy/hydro_io.h
index 59c49c0665f22eed3ef9b68de054e04533e80daa..f2247c234aa46ad881266417a6147db87f807afe 100644
--- a/src/hydro/PressureEntropy/hydro_io.h
+++ b/src/hydro/PressureEntropy/hydro_io.h
@@ -42,8 +42,9 @@
  * @param list The list of i/o properties to read.
  * @param num_fields The number of i/o fields to read.
  */
-void hydro_read_particles(struct part* parts, struct io_props* list,
-                          int* num_fields) {
+INLINE static void hydro_read_particles(struct part* parts,
+                                        struct io_props* list,
+                                        int* num_fields) {
 
   *num_fields = 8;
 
@@ -67,20 +68,21 @@ void hydro_read_particles(struct part* parts, struct io_props* list,
                                 UNIT_CONV_DENSITY, parts, rho);
 }
 
-void convert_u(const struct engine* e, const struct part* p,
-               const struct xpart* xp, float* ret) {
+INLINE static void convert_u(const struct engine* e, const struct part* p,
+                             const struct xpart* xp, float* ret) {
 
   ret[0] = hydro_get_comoving_internal_energy(p);
 }
 
-void convert_P(const struct engine* e, const struct part* p,
-               const struct xpart* xp, float* ret) {
+INLINE static void convert_P(const struct engine* e, const struct part* p,
+                             const struct xpart* xp, float* ret) {
 
   ret[0] = hydro_get_comoving_pressure(p);
 }
 
-void convert_part_pos(const struct engine* e, const struct part* p,
-                      const struct xpart* xp, double* ret) {
+INLINE static void convert_part_pos(const struct engine* e,
+                                    const struct part* p,
+                                    const struct xpart* xp, double* ret) {
 
   if (e->s->periodic) {
     ret[0] = box_wrap(p->x[0], 0.0, e->s->dim[0]);
@@ -93,8 +95,9 @@ void convert_part_pos(const struct engine* e, const struct part* p,
   }
 }
 
-void convert_part_vel(const struct engine* e, const struct part* p,
-                      const struct xpart* xp, float* ret) {
+INLINE static void convert_part_vel(const struct engine* e,
+                                    const struct part* p,
+                                    const struct xpart* xp, float* ret) {
 
   const int with_cosmology = (e->policy & engine_policy_cosmology);
   const struct cosmology* cosmo = e->cosmology;
@@ -127,8 +130,9 @@ void convert_part_vel(const struct engine* e, const struct part* p,
   ret[2] *= cosmo->a2_inv;
 }
 
-void convert_part_potential(const struct engine* e, const struct part* p,
-                            const struct xpart* xp, float* ret) {
+INLINE static void convert_part_potential(const struct engine* e,
+                                          const struct part* p,
+                                          const struct xpart* xp, float* ret) {
 
   if (p->gpart != NULL)
     ret[0] = gravity_get_comoving_potential(p->gpart);
@@ -143,8 +147,10 @@ void convert_part_potential(const struct engine* e, const struct part* p,
  * @param list The list of i/o properties to write.
  * @param num_fields The number of i/o fields to write.
  */
-void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
-                           struct io_props* list, int* num_fields) {
+INLINE static void hydro_write_particles(const struct part* parts,
+                                         const struct xpart* xparts,
+                                         struct io_props* list,
+                                         int* num_fields) {
 
   *num_fields = 11;
 
@@ -180,7 +186,7 @@ void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
  * @brief Writes the current model of SPH to the file
  * @param h_grpsph The HDF5 group in which to write
  */
-void hydro_write_flavour(hid_t h_grpsph) {
+INLINE static void hydro_write_flavour(hid_t h_grpsph) {
 
   /* Viscosity and thermal conduction */
   /* Nothing in this minimal model... */
@@ -201,6 +207,6 @@ void hydro_write_flavour(hid_t h_grpsph) {
  *
  * @return 1 if entropy is in 'internal energy', 0 otherwise.
  */
-int writeEntropyFlag(void) { return 0; }
+INLINE static int writeEntropyFlag(void) { return 0; }
 
 #endif /* SWIFT_PRESSURE_ENTROPY_HYDRO_IO_H */
diff --git a/src/hydro/Shadowswift/hydro_io.h b/src/hydro/Shadowswift/hydro_io.h
index 4a15ec845dd3913b10c1b6579290f2ec6bdb4e84..1f6bb86e62c6a3359d1242328775c6e4067ef8f2 100644
--- a/src/hydro/Shadowswift/hydro_io.h
+++ b/src/hydro/Shadowswift/hydro_io.h
@@ -32,8 +32,9 @@
  * @param list The list of i/o properties to read.
  * @param num_fields The number of i/o fields to read.
  */
-void hydro_read_particles(struct part* parts, struct io_props* list,
-                          int* num_fields) {
+INLINE static void hydro_read_particles(struct part* parts,
+                                        struct io_props* list,
+                                        int* num_fields) {
 
   *num_fields = 8;
 
@@ -64,8 +65,8 @@ void hydro_read_particles(struct part* parts, struct io_props* list,
  * @param p Particle.
  * @return Internal energy of the particle
  */
-void convert_u(const struct engine* e, const struct part* p,
-               const struct xpart* xp, float* ret) {
+INLINE static void convert_u(const struct engine* e, const struct part* p,
+                             const struct xpart* xp, float* ret) {
   ret[0] = hydro_get_internal_energy(p);
 }
 
@@ -76,8 +77,8 @@ void convert_u(const struct engine* e, const struct part* p,
  * @param p Particle.
  * @return Entropic function of the particle
  */
-void convert_A(const struct engine* e, const struct part* p,
-               const struct xpart* xp, float* ret) {
+INLINE static void convert_A(const struct engine* e, const struct part* p,
+                             const struct xpart* xp, float* ret) {
   ret[0] = hydro_get_entropy(p);
 }
 
@@ -88,8 +89,8 @@ void convert_A(const struct engine* e, const struct part* p,
  * @param p Particle.
  * @return Total energy of the particle
  */
-void convert_Etot(const struct engine* e, const struct part* p,
-                  const struct xpart* xp, float* ret) {
+INLINE static void convert_Etot(const struct engine* e, const struct part* p,
+                                const struct xpart* xp, float* ret) {
 #ifdef SHADOWFAX_TOTAL_ENERGY
   return p->conserved.energy;
 #else
@@ -107,8 +108,9 @@ void convert_Etot(const struct engine* e, const struct part* p,
 #endif
 }
 
-void convert_part_pos(const struct engine* e, const struct part* p,
-                      const struct xpart* xp, double* ret) {
+INLINE static void convert_part_pos(const struct engine* e,
+                                    const struct part* p,
+                                    const struct xpart* xp, double* ret) {
 
   if (e->s->periodic) {
     ret[0] = box_wrap(p->x[0], 0.0, e->s->dim[0]);
@@ -128,8 +130,10 @@ void convert_part_pos(const struct engine* e, const struct part* p,
  * @param list The list of i/o properties to write.
  * @param num_fields The number of i/o fields to write.
  */
-void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
-                           struct io_props* list, int* num_fields) {
+INLINE static void hydro_write_particles(const struct part* parts,
+                                         const struct xpart* xparts,
+                                         struct io_props* list,
+                                         int* num_fields) {
 
   *num_fields = 13;
 
@@ -168,7 +172,7 @@ void hydro_write_particles(const struct part* parts, const struct xpart* xparts,
  * @brief Writes the current model of SPH to the file
  * @param h_grpsph The HDF5 group in which to write
  */
-void hydro_write_flavour(hid_t h_grpsph) {
+INLINE static void hydro_write_flavour(hid_t h_grpsph) {
   /* Gradient information */
   io_write_attribute_s(h_grpsph, "Gradient reconstruction model",
                        HYDRO_GRADIENT_IMPLEMENTATION);
@@ -189,4 +193,4 @@ void hydro_write_flavour(hid_t h_grpsph) {
  *
  * @return 1 if entropy is in 'internal energy', 0 otherwise.
  */
-int writeEntropyFlag(void) { return 0; }
+INLINE static int writeEntropyFlag(void) { return 0; }
diff --git a/src/parallel_io.c b/src/parallel_io.c
index e09022516c22f870413538754b08ad8b911d27de..3c073ce712c4bdcc97dd70ea6e146d92f2f1b64d 100644
--- a/src/parallel_io.c
+++ b/src/parallel_io.c
@@ -49,6 +49,7 @@
 #include "io_properties.h"
 #include "kernel_hydro.h"
 #include "part.h"
+#include "part_type.h"
 #include "stars_io.h"
 #include "units.h"
 #include "xmf.h"
@@ -838,6 +839,7 @@ void prepare_file(struct engine* e, const char* baseName, long long N_total[6],
   const struct xpart* xparts = e->s->xparts;
   const struct gpart* gparts = e->s->gparts;
   const struct spart* sparts = e->s->sparts;
+  const struct swift_params* params = e->parameter_file;
   FILE* xmfFile = 0;
   int periodic = e->s->periodic;
   int numFiles = 1;
@@ -1017,10 +1019,19 @@ void prepare_file(struct engine* e, const char* baseName, long long N_total[6],
         error("Particle Type %d not yet supported. Aborting", ptype);
     }
 
-    /* Prepare everything */
-    for (int i = 0; i < num_fields; ++i)
-      prepareArray(e, h_grp, fileName, xmfFile, partTypeGroupName, list[i],
-                   N_total[ptype], snapshot_units);
+    /* Prepare everything that is not cancelled */
+    for (int i = 0; i < num_fields; ++i) {
+
+      /* Did the user cancel this field? */
+      char field[PARSER_MAX_LINE_SIZE];
+      sprintf(field, "SelectOutput:%s_%s", list[i].name,
+              part_type_names[ptype]);
+      int should_write = parser_get_opt_param_int(params, field, 1);
+
+      if (should_write)
+        prepareArray(e, h_grp, fileName, xmfFile, partTypeGroupName, list[i],
+                     N_total[ptype], snapshot_units);
+    }
 
     /* Close particle group */
     H5Gclose(h_grp);
@@ -1072,6 +1083,7 @@ void write_output_parallel(struct engine* e, const char* baseName,
   struct gpart* dmparts = NULL;
   const struct spart* sparts = e->s->sparts;
   const struct cooling_function_data* cooling = e->cooling_func;
+  const struct swift_params* params = e->parameter_file;
 
   /* Number of unassociated gparts */
   const size_t Ndm = Ntot > 0 ? Ntot - (Ngas + Nstars) : 0;
@@ -1261,11 +1273,20 @@ void write_output_parallel(struct engine* e, const char* baseName,
         error("Particle Type %d not yet supported. Aborting", ptype);
     }
 
-    /* Write everything */
-    for (int i = 0; i < num_fields; ++i)
-      writeArray(e, h_grp, fileName, partTypeGroupName, list[i], Nparticles,
-                 N_total[ptype], mpi_rank, offset[ptype], internal_units,
-                 snapshot_units);
+    /* Write everything that is not cancelled */
+    for (int i = 0; i < num_fields; ++i) {
+
+      /* Did the user cancel this field? */
+      char field[PARSER_MAX_LINE_SIZE];
+      sprintf(field, "SelectOutput:%s_%s", list[i].name,
+              part_type_names[ptype]);
+      int should_write = parser_get_opt_param_int(params, field, 1);
+
+      if (should_write)
+        writeArray(e, h_grp, fileName, partTypeGroupName, list[i], Nparticles,
+                   N_total[ptype], mpi_rank, offset[ptype], internal_units,
+                   snapshot_units);
+    }
 
     /* Free temporary array */
     if (dmparts) {
diff --git a/src/parser.c b/src/parser.c
index af9ef5fd6b7228dd4ad96640b59322175586862b..581f3a0cffd531b5c879e80dc6c256e6761666ee 100644
--- a/src/parser.c
+++ b/src/parser.c
@@ -57,12 +57,23 @@ static void find_duplicate_section(const struct swift_params *params,
 static int lineNumber = 0;
 
 /**
- * @brief Reads an input file and stores each parameter in a structure.
+ * @brief Initialize the parser structure.
  *
  * @param file_name Name of file to be read
  * @param params Structure to be populated from file
  */
+void parser_init(const char *file_name, struct swift_params *params) {
+  params->paramCount = 0;
+  params->sectionCount = 0;
+  strcpy(params->fileName, file_name);
+}
 
+/**
+ * @brief Reads an input file and stores each parameter in a structure.
+ *
+ * @param file_name Name of file to be read
+ * @param params Structure to be populated from file
+ */
 void parser_read_file(const char *file_name, struct swift_params *params) {
   /* Open file for reading */
   FILE *file = fopen(file_name, "r");
@@ -71,9 +82,7 @@ void parser_read_file(const char *file_name, struct swift_params *params) {
   char line[PARSER_MAX_LINE_SIZE];
 
   /* Initialise parameter count. */
-  params->paramCount = 0;
-  params->sectionCount = 0;
-  strcpy(params->fileName, file_name);
+  parser_init(file_name, params);
 
   /* Check if parameter file exits. */
   if (file == NULL) {
diff --git a/src/parser.h b/src/parser.h
index 58f4a53ea2114dbed29f49e2b5ccca69239b45e3..dc4520b9860798ac1d902a1f100346f22ec3b0e3 100644
--- a/src/parser.h
+++ b/src/parser.h
@@ -55,6 +55,7 @@ struct swift_params {
 };
 
 /* Public API. */
+void parser_init(const char *file_name, struct swift_params *params);
 void parser_read_file(const char *file_name, struct swift_params *params);
 void parser_print_params(const struct swift_params *params);
 void parser_write_params_to_file(const struct swift_params *params,
diff --git a/src/serial_io.c b/src/serial_io.c
index 9403caad7670b9af369f4b3598b8a05cf2d0d9e9..750ec23bfd486e5f4f6226d4790a0af5b306b873 100644
--- a/src/serial_io.c
+++ b/src/serial_io.c
@@ -49,6 +49,7 @@
 #include "io_properties.h"
 #include "kernel_hydro.h"
 #include "part.h"
+#include "part_type.h"
 #include "stars_io.h"
 #include "units.h"
 #include "xmf.h"
@@ -727,6 +728,7 @@ void write_output_serial(struct engine* e, const char* baseName,
   struct gpart* dmparts = NULL;
   const struct spart* sparts = e->s->sparts;
   const struct cooling_function_data* cooling = e->cooling_func;
+  const struct swift_params* params = e->parameter_file;
   FILE* xmfFile = 0;
 
   /* Number of unassociated gparts */
@@ -1005,11 +1007,20 @@ void write_output_serial(struct engine* e, const char* baseName,
             error("Particle Type %d not yet supported. Aborting", ptype);
         }
 
-        /* Write everything */
-        for (int i = 0; i < num_fields; ++i)
-          writeArray(e, h_grp, fileName, xmfFile, partTypeGroupName, list[i],
-                     Nparticles, N_total[ptype], mpi_rank, offset[ptype],
-                     internal_units, snapshot_units);
+        /* Write everything that is not cancelled */
+        for (int i = 0; i < num_fields; ++i) {
+
+          /* Did the user cancel this field? */
+          char field[PARSER_MAX_LINE_SIZE];
+          sprintf(field, "SelectOutput:%s_%s", list[i].name,
+                  part_type_names[ptype]);
+          int should_write = parser_get_opt_param_int(params, field, 1);
+
+          if (should_write)
+            writeArray(e, h_grp, fileName, xmfFile, partTypeGroupName, list[i],
+                       Nparticles, N_total[ptype], mpi_rank, offset[ptype],
+                       internal_units, snapshot_units);
+        }
 
         /* Free temporary array */
         if (dmparts) {
diff --git a/src/single_io.c b/src/single_io.c
index 866d830ac59f2f7f6d7e84212f80e715ed066b6c..47599eaaa39de21a3df6ad810b7dd71b4807c995 100644
--- a/src/single_io.c
+++ b/src/single_io.c
@@ -48,6 +48,7 @@
 #include "io_properties.h"
 #include "kernel_hydro.h"
 #include "part.h"
+#include "part_type.h"
 #include "stars_io.h"
 #include "units.h"
 #include "xmf.h"
@@ -594,6 +595,7 @@ void write_output_single(struct engine* e, const char* baseName,
   struct gpart* dmparts = NULL;
   const struct spart* sparts = e->s->sparts;
   const struct cooling_function_data* cooling = e->cooling_func;
+  const struct swift_params* params = e->parameter_file;
 
   /* Number of unassociated gparts */
   const size_t Ndm = Ntot > 0 ? Ntot - (Ngas + Nstars) : 0;
@@ -824,10 +826,19 @@ void write_output_single(struct engine* e, const char* baseName,
         error("Particle Type %d not yet supported. Aborting", ptype);
     }
 
-    /* Write everything */
-    for (int i = 0; i < num_fields; ++i)
-      writeArray(e, h_grp, fileName, xmfFile, partTypeGroupName, list[i], N,
-                 internal_units, snapshot_units);
+    /* Write everything that is not cancelled */
+    for (int i = 0; i < num_fields; ++i) {
+
+      /* Did the user cancel this field? */
+      char field[PARSER_MAX_LINE_SIZE];
+      sprintf(field, "SelectOutput:%s_%s", list[i].name,
+              part_type_names[ptype]);
+      int should_write = parser_get_opt_param_int(params, field, 1);
+
+      if (should_write)
+        writeArray(e, h_grp, fileName, xmfFile, partTypeGroupName, list[i], N,
+                   internal_units, snapshot_units);
+    }
 
     /* Free temporary array */
     if (dmparts) {
diff --git a/src/stars/Default/star_io.h b/src/stars/Default/star_io.h
index c3dc31096383533e1e15fa65615d2c9aac0f43e3..7ad29f0a935c002b1337c2a75d6f987c05c9bb43 100644
--- a/src/stars/Default/star_io.h
+++ b/src/stars/Default/star_io.h
@@ -28,8 +28,8 @@
  * @param list The list of i/o properties to read.
  * @param num_fields The number of i/o fields to read.
  */
-void star_read_particles(struct spart* sparts, struct io_props* list,
-                         int* num_fields) {
+INLINE static void star_read_particles(struct spart* sparts,
+                                       struct io_props* list, int* num_fields) {
 
   /* Say how much we want to read */
   *num_fields = 4;
@@ -52,8 +52,9 @@ void star_read_particles(struct spart* sparts, struct io_props* list,
  * @param list The list of i/o properties to write.
  * @param num_fields The number of i/o fields to write.
  */
-void star_write_particles(const struct spart* sparts, struct io_props* list,
-                          int* num_fields) {
+INLINE static void star_write_particles(const struct spart* sparts,
+                                        struct io_props* list,
+                                        int* num_fields) {
 
   /* Say how much we want to read */
   *num_fields = 4;
diff --git a/src/swift.h b/src/swift.h
index 7691720942f32d29d3269ccdd7adbe8db32280bf..17f7c5c7e4e56ed1650d6cdd078060244c67b643 100644
--- a/src/swift.h
+++ b/src/swift.h
@@ -29,8 +29,10 @@
 #include "cell.h"
 #include "chemistry.h"
 #include "clocks.h"
+#include "common_io.h"
 #include "const.h"
 #include "cooling.h"
+#include "cooling_struct.h"
 #include "cosmology.h"
 #include "cycle.h"
 #include "debug.h"
diff --git a/src/version.c b/src/version.c
index a0a2ae42a5ce1a24ffce43f7ae3f59e1d6ed4168..fc0a8653a4c31e3bb57d5894f2e812fcaa6bcc2d 100644
--- a/src/version.c
+++ b/src/version.c
@@ -146,7 +146,7 @@ const char *configuration_options(void) {
   static int initialised = 0;
   static const char *config = SWIFT_CONFIG_FLAGS;
   if (!initialised) {
-    snprintf(buf, 1024, "'%s'", config);
+    snprintf(buf, 1024, "'%.1022s'", config);
     initialised = 1;
   }
   return buf;
@@ -162,7 +162,7 @@ const char *compilation_cflags(void) {
   static int initialised = 0;
   static const char *cflags = SWIFT_CFLAGS;
   if (!initialised) {
-    snprintf(buf, 1024, "'%s'", cflags);
+    snprintf(buf, 1024, "'%.1022s'", cflags);
     initialised = 1;
   }
   return buf;
@@ -272,12 +272,12 @@ const char *mpi_version(void) {
 #else
   /* Use autoconf guessed value. */
   static char lib_version[60] = {0};
-  snprintf(lib_version, 60, "%s", SWIFT_MPI_LIBRARY);
+  snprintf(lib_version, 60, "%.60s", SWIFT_MPI_LIBRARY);
 #endif
 
   /* Numeric version. */
   MPI_Get_version(&std_version, &std_subversion);
-  snprintf(version, 80, "%s (MPI std v%i.%i)", lib_version, std_version,
+  snprintf(version, 80, "%.60s (MPI std v%i.%i)", lib_version, std_version,
            std_subversion);
 #else
   sprintf(version, "Code was not compiled with MPI support");
@@ -345,7 +345,7 @@ const char *libgsl_version(void) {
 
   static char version[256] = {0};
 #if defined(HAVE_LIBGSL)
-  sprintf(version, "%s", gsl_version);
+  sprintf(version, "%.255s", gsl_version);
 #else
   sprintf(version, "Unknown version");
 #endif
diff --git a/tests/Makefile.am b/tests/Makefile.am
index 564076ed501409083db017de4c79c17f3b2a7c15..e8f2a4f0ae39254c28e981fef5e36866e56395d8 100644
--- a/tests/Makefile.am
+++ b/tests/Makefile.am
@@ -27,7 +27,7 @@ TESTS = testGreetings testMaths testReading.sh testSingle testKernel testSymmetr
         testMatrixInversion testThreadpool testDump testLogger testInteractions.sh \
         testVoronoi1D testVoronoi2D testVoronoi3D testGravityDerivatives \
 	testPeriodicBC.sh testPeriodicBCPerturbed.sh testPotentialSelf \
-	testPotentialPair testEOS testUtilities
+	testPotentialPair testEOS testUtilities testSelectOutput.sh
 
 # List of test programs to compile
 check_PROGRAMS = testGreetings testReading testSingle testTimeIntegration \
@@ -37,7 +37,8 @@ check_PROGRAMS = testGreetings testReading testSingle testTimeIntegration \
                  testAdiabaticIndex testRiemannExact testRiemannTRRS \
                  testRiemannHLLC testMatrixInversion testDump testLogger \
 		 testVoronoi1D testVoronoi2D testVoronoi3D testPeriodicBC \
-		 testGravityDerivatives testPotentialSelf testPotentialPair testEOS testUtilities
+		 testGravityDerivatives testPotentialSelf testPotentialPair testEOS testUtilities \
+		 testSelectOutput
 
 # Rebuild tests when SWIFT is updated.
 $(check_PROGRAMS): ../src/.libs/libswiftsim.a
@@ -49,6 +50,8 @@ testMaths_SOURCES = testMaths.c
 
 testReading_SOURCES = testReading.c
 
+testSelectOutput_SOURCES = testSelectOutput.c
+
 testSymmetry_SOURCES = testSymmetry.c
 
 # Added because of issues using memcmp on clang 4.x
@@ -120,4 +123,4 @@ EXTRA_DIST = testReading.sh makeInput.py testActivePair.sh \
              tolerance_27_normal.dat tolerance_27_perturbed.dat tolerance_27_perturbed_h.dat tolerance_27_perturbed_h2.dat \
 	     tolerance_testInteractions.dat tolerance_pair_active.dat tolerance_pair_force_active.dat \
 	     fft_params.yml tolerance_periodic_BC_normal.dat tolerance_periodic_BC_perturbed.dat \
-	     testEOS.sh testEOS_plot.sh
+	     testEOS.sh testEOS_plot.sh testSelectOutput.sh selectOutput.yml
diff --git a/tests/selectOutput.yml b/tests/selectOutput.yml
new file mode 100644
index 0000000000000000000000000000000000000000..1778935146b19992e25efcb320d8cc523c6472a5
--- /dev/null
+++ b/tests/selectOutput.yml
@@ -0,0 +1,12 @@
+SelectOutput:
+  # Particle Type 0
+  Coordinates_Gas: 1   # check if written when specified
+  Velocities_Gas:  0   # check if not written when specified
+  Masses_Gas:     -5   # check warning if not 0 or 1 and if written
+  Pot_Gas:         1   # check warning if wrong name
+  # Density_Gas:   1   # check if written when not specified
+
+# Parameters for the hydrodynamics scheme
+SPH:
+  resolution_eta:        1.2348   # Target smoothing length in units of the mean inter-particle separation (1.2348 == 48Ngbs with the cubic spline kernel).
+  CFL_condition:         0.1      # Courant-Friedrich-Levy condition for time integration.
diff --git a/tests/test125cells.c b/tests/test125cells.c
index 41023aec16d1a3619d58caf8fbb01f08285f0661..ddce3176d463f2ef754bf82b364858085e317e4b 100644
--- a/tests/test125cells.c
+++ b/tests/test125cells.c
@@ -759,7 +759,7 @@ int main(int argc, char *argv[]) {
 
     /* Dump if necessary */
     if (n == 0) {
-      sprintf(outputFileName, "swift_dopair_125_%s.dat",
+      sprintf(outputFileName, "swift_dopair_125_%.150s.dat",
               outputFileNameExtension);
       dump_particle_fields(outputFileName, main_cell, solution, 0);
     }
@@ -876,7 +876,8 @@ int main(int argc, char *argv[]) {
   /* Output timing */
   message("Brute force calculation took : %15lli ticks.", toc - tic);
 
-  sprintf(outputFileName, "brute_force_125_%s.dat", outputFileNameExtension);
+  sprintf(outputFileName, "brute_force_125_%.150s.dat",
+          outputFileNameExtension);
   dump_particle_fields(outputFileName, main_cell, solution, 0);
 
   /* Clean things to make the sanitizer happy ... */
diff --git a/tests/test27cells.c b/tests/test27cells.c
index ccaaad8e8644ae72973c9b93eb0af44c77061113..ada1b782cfff3866bf26937391007947e9c9a175 100644
--- a/tests/test27cells.c
+++ b/tests/test27cells.c
@@ -546,7 +546,7 @@ int main(int argc, char *argv[]) {
 
     /* Dump if necessary */
     if (i % 50 == 0) {
-      sprintf(outputFileName, "swift_dopair_27_%s.dat",
+      sprintf(outputFileName, "swift_dopair_27_%.150s.dat",
               outputFileNameExtension);
       dump_particle_fields(outputFileName, main_cell, cells);
     }
@@ -589,7 +589,7 @@ int main(int argc, char *argv[]) {
   end_calculation(main_cell, &cosmo);
 
   /* Dump */
-  sprintf(outputFileName, "brute_force_27_%s.dat", outputFileNameExtension);
+  sprintf(outputFileName, "brute_force_27_%.150s.dat", outputFileNameExtension);
   dump_particle_fields(outputFileName, main_cell, cells);
 
   /* Output timing */
diff --git a/tests/testActivePair.c b/tests/testActivePair.c
index 4a70e341d0160193da595274f47e943c76390c15..6889a18887894af0a9434f786df21dbf842e87e5 100644
--- a/tests/testActivePair.c
+++ b/tests/testActivePair.c
@@ -578,8 +578,9 @@ int main(int argc, char *argv[]) {
   runner->e = &engine;
 
   /* Create output file names. */
-  sprintf(swiftOutputFileName, "swift_dopair_%s.dat", outputFileNameExtension);
-  sprintf(bruteForceOutputFileName, "brute_force_pair_%s.dat",
+  sprintf(swiftOutputFileName, "swift_dopair_%.150s.dat",
+          outputFileNameExtension);
+  sprintf(bruteForceOutputFileName, "brute_force_pair_%.150s.dat",
           outputFileNameExtension);
 
   /* Delete files if they already exist. */
@@ -632,9 +633,9 @@ int main(int argc, char *argv[]) {
   finalise = &end_calculation_force;
 
   /* Create new output file names. */
-  sprintf(swiftOutputFileName, "swift_dopair2_force_%s.dat",
+  sprintf(swiftOutputFileName, "swift_dopair2_force_%.150s.dat",
           outputFileNameExtension);
-  sprintf(bruteForceOutputFileName, "brute_force_dopair2_%s.dat",
+  sprintf(bruteForceOutputFileName, "brute_force_dopair2_%.150s.dat",
           outputFileNameExtension);
 
   /* Delete files if they already exist. */
diff --git a/tests/testPeriodicBC.c b/tests/testPeriodicBC.c
index 3e4cb5d3816ad706cd5f083912da7f453a2a0e5a..ffaa3bda0ccb62cd44169e228086267d2399c31f 100644
--- a/tests/testPeriodicBC.c
+++ b/tests/testPeriodicBC.c
@@ -512,9 +512,9 @@ int main(int argc, char *argv[]) {
   }
 
   /* Create output file names. */
-  sprintf(swiftOutputFileName, "swift_periodic_BC_%s.dat",
+  sprintf(swiftOutputFileName, "swift_periodic_BC_%.150s.dat",
           outputFileNameExtension);
-  sprintf(bruteForceOutputFileName, "brute_force_periodic_BC_%s.dat",
+  sprintf(bruteForceOutputFileName, "brute_force_periodic_BC_%.150s.dat",
           outputFileNameExtension);
 
   /* Delete files if they already exist. */
diff --git a/tests/testSelectOutput.c b/tests/testSelectOutput.c
new file mode 100644
index 0000000000000000000000000000000000000000..3bedddd03784bafddd4599c124929a6d88fbd9b0
--- /dev/null
+++ b/tests/testSelectOutput.c
@@ -0,0 +1,162 @@
+/*******************************************************************************
+ * This file is part of SWIFT.
+ * Copyright (C) 2015 Matthieu Schaller (matthieu.schaller@durham.ac.uk).
+ *
+ * This program 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 3 of the License, or
+ * (at your option) any later version.
+ *
+ * This program 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 General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public License
+ * along with this program.  If not, see <http://www.gnu.org/licenses/>.
+ *
+ ******************************************************************************/
+
+/* Some standard headers. */
+#include <stdlib.h>
+
+/* Includes. */
+#include "swift.h"
+
+void select_output_engine_init(struct engine *e, struct space *s,
+                               struct cosmology *cosmo,
+                               struct swift_params *params,
+                               struct cooling_function_data *cooling,
+                               struct hydro_props *hydro_properties) {
+  /* set structures */
+  e->s = s;
+  e->cooling_func = cooling;
+  e->parameter_file = params;
+  e->cosmology = cosmo;
+  e->policy = engine_policy_hydro;
+  e->hydro_properties = hydro_properties;
+
+  /* initialization of threadpool */
+  threadpool_init(&e->threadpool, 1);
+
+  /* set parameters */
+  e->verbose = 1;
+  e->time = 0;
+  e->snapshot_output_count = 0;
+  e->snapshot_compression = 0;
+  e->snapshot_label_delta = 1;
+};
+
+void select_output_space_init(struct space *s, double *dim, int periodic,
+                              size_t Ngas, size_t Nspart, size_t Ngpart,
+                              struct part *parts, struct spart *sparts,
+                              struct gpart *gparts) {
+  s->periodic = periodic;
+  for (int i = 0; i < 3; i++) {
+    s->dim[i] = dim[i];
+  }
+
+  /* init space particles */
+  s->nr_parts = Ngas;
+  s->nr_sparts = Nspart;
+  s->nr_gparts = Ngpart;
+
+  s->parts = parts;
+  s->gparts = gparts;
+  s->sparts = sparts;
+
+  /* Allocate the extra parts array for the gas particles. */
+  if (posix_memalign((void **)&s->xparts, xpart_align,
+                     Ngas * sizeof(struct xpart)) != 0)
+    error("Failed to allocate xparts.");
+  bzero(s->xparts, Ngas * sizeof(struct xpart));
+};
+
+void select_output_space_clean(struct space *s) { free(s->xparts); };
+
+void select_output_engine_clean(struct engine *e) {
+  threadpool_clean(&e->threadpool);
+}
+
+int main(int argc, char *argv[]) {
+
+  /* Initialize CPU frequency, this also starts time. */
+  unsigned long long cpufreq = 0;
+  clocks_set_cpufreq(cpufreq);
+
+  char *base_name = "testSelectOutput";
+  size_t Ngas = 0, Ngpart = 0, Nspart = 0;
+  int periodic = -1;
+  int flag_entropy_ICs = -1;
+  double dim[3];
+  struct part *parts = NULL;
+  struct gpart *gparts = NULL;
+  struct spart *sparts = NULL;
+
+  /* parse parameters */
+  message("Reading parameters.");
+  struct swift_params param_file;
+  char *input_file = "selectOutput.yml";
+  parser_read_file(input_file, &param_file);
+
+  /* Default unit system */
+  message("Initialization of the unit system.");
+  struct unit_system us;
+  units_init_cgs(&us);
+
+  /* Default physical constants */
+  message("Initialization of the physical constants.");
+  struct phys_const prog_const;
+  phys_const_init(&us, &param_file, &prog_const);
+
+  /* Read data */
+  message("Reading initial conditions.");
+  read_ic_single("input.hdf5", &us, dim, &parts, &gparts, &sparts, &Ngas,
+                 &Ngpart, &Nspart, &periodic, &flag_entropy_ICs, 1, 0, 0, 0, 1.,
+                 1, 0);
+
+  /* pseudo initialization of the space */
+  message("Initialization of the space.");
+  struct space s;
+  select_output_space_init(&s, dim, periodic, Ngas, Nspart, Ngpart, parts,
+                           sparts, gparts);
+
+  /* initialization of cosmology */
+  message("Initialization of the cosmology.");
+  struct cosmology cosmo;
+  cosmology_init_no_cosmo(&cosmo);
+
+  /* pseudo initialization of cooling */
+  message("Initialization of the cooling.");
+  struct cooling_function_data cooling;
+
+  /* pseudo initialization of hydro */
+  message("Initialization of the hydro.");
+  struct hydro_props hydro_properties;
+  hydro_props_init(&hydro_properties, &prog_const, &us, &param_file);
+
+  /* pseudo initialization of the engine */
+  message("Initialization of the engine.");
+  struct engine e;
+  select_output_engine_init(&e, &s, &cosmo, &param_file, &cooling,
+                            &hydro_properties);
+
+  /* check output selection */
+  message("Checking output parameters.");
+  long long N_total[swift_type_count] = {Ngas, Ngpart, 0, 0, Nspart, 0};
+  io_check_output_fields(&param_file, N_total);
+
+  /* write output file */
+  message("Writing output.");
+  write_output_single(&e, base_name, &us, &us);
+
+  /* Clean-up */
+  message("Cleaning memory.");
+  select_output_engine_clean(&e);
+  select_output_space_clean(&s);
+  cosmology_clean(&cosmo);
+  free(parts);
+  free(gparts);
+
+  return 0;
+}
diff --git a/tests/testSelectOutput.py b/tests/testSelectOutput.py
new file mode 100644
index 0000000000000000000000000000000000000000..aec7f4671fb2768acde768fd9929168559ebb3cb
--- /dev/null
+++ b/tests/testSelectOutput.py
@@ -0,0 +1,54 @@
+###############################################################################
+ # This file is part of SWIFT.
+ # Copyright (c) 2015 Bert Vandenbroucke (bert.vandenbroucke@ugent.be)
+ #                    Matthieu Schaller (matthieu.schaller@durham.ac.uk)
+ # 
+ # This program 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 3 of the License, or
+ # (at your option) any later version.
+ # 
+ # This program 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 General Public License for more details.
+ # 
+ # You should have received a copy of the GNU Lesser General Public License
+ # along with this program.  If not, see <http://www.gnu.org/licenses/>.
+ # 
+ ##############################################################################
+
+# Check the output done with swift
+
+import h5py
+
+filename = "testSelectOutput_0000.hdf5"
+log_filename = "select_output.log"
+
+# Read the simulation data
+sim = h5py.File(filename, "r")
+part0 = sim["/PartType0"]
+
+# check presence / absence fields
+if "Velocities" in part0:
+    raise Exception("`Velocities` present in HDF5 but should not be written")
+
+if "Coordinates" not in part0:
+    raise Exception("`Coordinates` not present in HDF5 but should be written")
+
+if "Masses" not in part0:
+    raise Exception("`Masses` not present in HDF5 but should be written")
+
+if "Density" not in part0:
+    raise Exception("`Density` not present in HDF5 but should be written")
+
+
+# check error detection
+with open(log_filename, "r") as f:
+    data = f.read()
+
+if "SelectOutput:Masses_Gas" not in data:
+    raise Exception("Input error in `SelectOuput:Masses_Gas` not detected")
+
+if "SelectOutput:Pot_Gas" not in data:
+    raise Exception("Parameter name error not detected for `SelectOutput:Pot_Gas`")
diff --git a/tests/testSelectOutput.sh.in b/tests/testSelectOutput.sh.in
new file mode 100644
index 0000000000000000000000000000000000000000..85fd999643f82fd10d96013ad360a75a441a9e1a
--- /dev/null
+++ b/tests/testSelectOutput.sh.in
@@ -0,0 +1,14 @@
+#!/bin/bash
+
+echo "Creating initial conditions"
+python @srcdir@/makeInput.py
+
+echo "Generating output"
+./testSelectOutput 2>&1 | tee select_output.log
+
+echo "Checking output"
+python @srcdir@/testSelectOutput.py
+
+rm -f testSelectOutput_0000.hdf5 testSelectOutput.xmf select_output.log
+
+echo "Test passed"