diff --git a/examples/CoolingHalo/cooling_halo.yml b/examples/CoolingHalo/cooling_halo.yml
index e4e8750f46bf7f46c692c4ccd3afe8453e0f606a..5cc77d2e01a8a50f9df79e1196a2ca356cfa36f4 100644
--- a/examples/CoolingHalo/cooling_halo.yml
+++ b/examples/CoolingHalo/cooling_halo.yml
@@ -31,15 +31,11 @@ SPH:
 # Parameters related to the initial conditions
 InitialConditions:
   file_name:  CoolingHalo.hdf5       # The file to read
-  shift_x:    0.                  # A shift to apply to all particles read from the ICs (in internal units).
-  shift_y:    0.
-  shift_z:    0.
+  shift:      [0.,0.,0.]             # A shift to apply to all particles read from the ICs (in internal units).
  
 # External potential parameters
 IsothermalPotential:
-  position_x:      0.     # location of centre of isothermal potential in internal units
-  position_y:      0.
-  position_z:      0.	
+  position:        [0.,0.,0.] # location of centre of isothermal potential in internal units
   vrot:            200.     # rotation speed of isothermal potential in internal units
   timestep_mult:   0.03     # controls time step
   epsilon:         0.1    #softening for the isothermal potential
diff --git a/examples/CoolingHaloWithSpin/cooling_halo.yml b/examples/CoolingHaloWithSpin/cooling_halo.yml
index 1d548b2fb7e531b712548eb5834b2c4de575f941..ddbd5d0b8498b52092fdb1f956eec6255249a39d 100644
--- a/examples/CoolingHaloWithSpin/cooling_halo.yml
+++ b/examples/CoolingHaloWithSpin/cooling_halo.yml
@@ -34,9 +34,7 @@ InitialConditions:
  
 # External potential parameters
 IsothermalPotential:
-  position_x:      0.     # Location of centre of isothermal potential in internal units
-  position_y:      0.
-  position_z:      0.	
+  position:        [0.,0.,0.] # Location of centre of isothermal potential in internal units
   vrot:            200.   # Rotation speed of isothermal potential in internal units
   timestep_mult:   0.03   # Controls time step
   epsilon:         1.0    # Softening for the isothermal potential
diff --git a/examples/ExternalPointMass/externalPointMass.yml b/examples/ExternalPointMass/externalPointMass.yml
index 6e2315cc4dd5c2c354375d39e03f13a1553ed08b..79dcd6458489958622ed7d6167d7f5f214ebd8b9 100644
--- a/examples/ExternalPointMass/externalPointMass.yml
+++ b/examples/ExternalPointMass/externalPointMass.yml
@@ -31,15 +31,11 @@ SPH:
 # Parameters related to the initial conditions
 InitialConditions:
   file_name:  PointMass.hdf5        # The file to read
-  shift_x:    50.                   # A shift to apply to all particles read from the ICs (in internal units).
-  shift_y:    50.
-  shift_z:    50.
+  shift:      [50.,50.,50.]         # A shift to apply to all particles read from the ICs (in internal units).
 
 # External potential parameters
 PointMassPotential:
-  position_x:      50.     # location of external point mass in internal units
-  position_y:      50.
-  position_z:      50.	
+  position:        [50.,50.,50.]    # location of external point mass in internal units
   mass:            1e10     # mass of external point mass in internal units
   timestep_mult:   0.03     # controls time step
 
diff --git a/examples/HydrostaticHalo/hydrostatic.yml b/examples/HydrostaticHalo/hydrostatic.yml
index 38a5df2866d39523e6e4b9f6c1c97a7537bff5a4..93dd80c02df283deb6a074485a3adb9023f9cfc8 100644
--- a/examples/HydrostaticHalo/hydrostatic.yml
+++ b/examples/HydrostaticHalo/hydrostatic.yml
@@ -34,9 +34,7 @@ InitialConditions:
  
 # External potential parameters
 IsothermalPotential:
-  position_x:      0.     # location of centre of isothermal potential in internal units
-  position_y:      0.
-  position_z:      0.	
+  position:        [0.,0.,0.] # location of centre of isothermal potential in internal units
   vrot:            200.     # rotation speed of isothermal potential in internal units
   epsilon:         1.0
   timestep_mult:   0.03     # controls time step
diff --git a/examples/IsothermalPotential/isothermal.yml b/examples/IsothermalPotential/isothermal.yml
index bd59ddc5117f71c478ce5eeda0a08817eefe8624..fd30895ceb617c676c497488b2ad379d63f5c03d 100644
--- a/examples/IsothermalPotential/isothermal.yml
+++ b/examples/IsothermalPotential/isothermal.yml
@@ -26,15 +26,11 @@ Snapshots:
 # Parameters related to the initial conditions
 InitialConditions:
   file_name:  Isothermal.hdf5       # The file to read
-  shift_x:    200.                  # Shift all particles to be in the potential
-  shift_y:    200.
-  shift_z:    200.
+  shift:      [200.,200.,200.]      # Shift all particles to be in the potential
  
 # External potential parameters
 IsothermalPotential:
-  position_x:      0.       # location of centre of isothermal potential in internal units
-  position_y:      0.
-  position_z:      0.
-  vrot:            200.     # rotation speed of isothermal potential in internal units
-  timestep_mult:   0.01     # controls time step
-  epsilon:         0.       # No softening at the centre of the halo
+  position:        [0.,0.,0.] # location of centre of isothermal potential in internal units
+  vrot:            200.       # rotation speed of isothermal potential in internal units
+  timestep_mult:   0.01       # controls time step
+  epsilon:         0.         # No softening at the centre of the halo
diff --git a/examples/KeplerianRing/keplerian_ring.yml b/examples/KeplerianRing/keplerian_ring.yml
index 421d1e89255195cd05a66975d838c59a4ad77c72..72a45eabdf9d217a35b440146722ca969eb50c60 100644
--- a/examples/KeplerianRing/keplerian_ring.yml
+++ b/examples/KeplerianRing/keplerian_ring.yml
@@ -32,15 +32,11 @@ SPH:
 # Parameters related to the initial conditions
 InitialConditions:
   file_name:  initial_conditions.hdf5        # The file to read
-  shift_x:    0.                   # A shift to apply to all particles read from the ICs (in internal units).
-  shift_y:    0.
-  shift_z:    0.
+  shift:      [0.,0.,0.]                     # A shift to apply to all particles read from the ICs (in internal units).
 
 # External potential parameters
 PointMassPotential:
-  position_x:      5.     # location of external point mass in internal units
-  position_y:      5.     # here just take this to be the centre of the ring
-  position_z:      5.	
-  mass:            1.498351e7     # mass of external point mass in internal units
-  timestep_mult:   0.03     # controls time step
+  position:        [5.,5.,5.]  # location of external point mass in internal units
+  mass:            1.498351e7  # mass of external point mass in internal units
+  timestep_mult:   0.03        # controls time step
   softening:       0.05
diff --git a/examples/KeplerianRing/write_gadget.py b/examples/KeplerianRing/write_gadget.py
index d2519cdaf4d78d9c7327419283072b8001e70fcb..0f745782ef81c7500e8b965938eb5eea33ede288 100644
--- a/examples/KeplerianRing/write_gadget.py
+++ b/examples/KeplerianRing/write_gadget.py
@@ -300,7 +300,7 @@ def write_block(f, part_type, pos, vel, ids, mass, int_energy, smoothing, other=
     else:
         data = default_data
 
-    particles = f.create_group(f"PartType{part_type}")
+    particles = f.create_group("PartType" + str(part_type))
 
     for name, value in data.items():
         particles.create_dataset(name, data=value)
diff --git a/examples/MultiTypes/multiTypes.yml b/examples/MultiTypes/multiTypes.yml
index 5872634cd620e0b7e2b7edb3a38e363c34229b9a..04647f0f00e69f5baf2560aca0feeb14a26cc50a 100644
--- a/examples/MultiTypes/multiTypes.yml
+++ b/examples/MultiTypes/multiTypes.yml
@@ -35,8 +35,6 @@ InitialConditions:
 
 # External potential parameters
 PointMassPotential:
-  position_x:      50.     # location of external point mass in internal units
-  position_y:      50.
-  position_z:      50.	
+  position:        [50.,50.,50.] # location of external point mass in internal units
   mass:            1e10     # mass of external point mass in internal units
   timestep_mult:   1e-2
diff --git a/examples/SodShock_3D/sodShock.yml b/examples/SodShock_3D/sodShock.yml
index 6042c8090d00fef5467a7fed3d6f5a104c626f43..d4234f2240adb7c9c956701476c884a9a55550e0 100644
--- a/examples/SodShock_3D/sodShock.yml
+++ b/examples/SodShock_3D/sodShock.yml
@@ -32,4 +32,5 @@ SPH:
 # Parameters related to the initial conditions
 InitialConditions:
   file_name:  ./sodShock.hdf5       # The file to read
+  shift: [1.,1.,0.]
 
diff --git a/examples/parameter_example.yml b/examples/parameter_example.yml
index 8b77cb9db84949eb4fc12c2c6bccca67f92cc899..78fd80a334b1570d1f091a8fab80698614207741 100644
--- a/examples/parameter_example.yml
+++ b/examples/parameter_example.yml
@@ -95,9 +95,7 @@ InitialConditions:
   cleanup_h_factors:           0    # (Optional) Clean up the h-factors used in the ICs (e.g. in Gadget files).
   cleanup_smoothing_lengths:   0    # (Optional) Clean the values of the smoothing lengths that are read in to remove stupid values. Set to 1 to activate.
   smoothing_length_scaling:    1.   # (Optional) A scaling factor to apply to all smoothing lengths in the ICs.
-  shift_x:    0.                    # (Optional) A shift to apply to all particles read from the ICs (in internal units).
-  shift_y:    0.
-  shift_z:    0.
+  shift:      [0.0,0.0,0.0]         # (Optional) A shift to apply to all particles read from the ICs (in internal units).
   replicate:  2                     # (Optional) Replicate all particles along each axis a given integer number of times. Default 1.
 
 # Parameters controlling restarts
@@ -145,17 +143,14 @@ EoS:
 
 # Point mass external potentials
 PointMassPotential:
-  position_x:      50.      # location of external point mass (internal units)
-  position_y:      50.
-  position_z:      50.
-  mass:            1e10     # mass of external point mass (internal units)
-  timestep_mult:   0.03     # Dimensionless pre-factor for the time-step condition
+  position:        [50.,50.0,50.]      # location of external point mass (internal units)
+  mass:            1e10                # mass of external point mass (internal units)
+  timestep_mult:   0.03                # Dimensionless pre-factor for the time-step condition
+  softening:       0.05                # For point-mass-softened option
 
 # Isothermal potential parameters
 IsothermalPotential:
-  position_x:      100.     # Location of centre of isothermal potential with respect to centre of the box (internal units)
-  position_y:      100.
-  position_z:      100.
+  position:        [100.,100.,100.]    # Location of centre of isothermal potential with respect to centre of the box (internal units)
   vrot:            200.     # Rotation speed of isothermal potential (internal units)
   timestep_mult:   0.03     # Dimensionless pre-factor for the time-step condition
   epsilon:         0.1      # Softening size (internal units)
diff --git a/src/potential/isothermal/potential.h b/src/potential/isothermal/potential.h
index 4267c9becc7251ffe276b69f078e374504c22aab..b5f8d7c39738bfe1895c73e6e59ae1279c0f74fa 100644
--- a/src/potential/isothermal/potential.h
+++ b/src/potential/isothermal/potential.h
@@ -42,7 +42,7 @@
 struct external_potential {
 
   /*! Position of the centre of potential */
-  double x, y, z;
+  double x[3];
 
   /*! Rotation velocity */
   double vrot;
@@ -73,9 +73,9 @@ __attribute__((always_inline)) INLINE static float external_gravity_timestep(
     const struct phys_const* restrict phys_const,
     const struct gpart* restrict g) {
 
-  const float dx = g->x[0] - potential->x;
-  const float dy = g->x[1] - potential->y;
-  const float dz = g->x[2] - potential->z;
+  const float dx = g->x[0] - potential->x[0];
+  const float dy = g->x[1] - potential->x[1];
+  const float dz = g->x[2] - potential->x[2];
 
   const float r2_plus_epsilon2_inv =
       1.f / (dx * dx + dy * dy + dz * dz + potential->epsilon2);
@@ -115,9 +115,9 @@ __attribute__((always_inline)) INLINE static void external_gravity_acceleration(
     double time, const struct external_potential* potential,
     const struct phys_const* const phys_const, struct gpart* g) {
 
-  const float dx = g->x[0] - potential->x;
-  const float dy = g->x[1] - potential->y;
-  const float dz = g->x[2] - potential->z;
+  const float dx = g->x[0] - potential->x[0];
+  const float dy = g->x[1] - potential->x[1];
+  const float dz = g->x[2] - potential->x[2];
   const float r2_plus_epsilon2_inv =
       1.f / (dx * dx + dy * dy + dz * dz + potential->epsilon2);
 
@@ -144,9 +144,9 @@ external_gravity_get_potential_energy(
     double time, const struct external_potential* potential,
     const struct phys_const* const phys_const, const struct gpart* g) {
 
-  const float dx = g->x[0] - potential->x;
-  const float dy = g->x[1] - potential->y;
-  const float dz = g->x[2] - potential->z;
+  const float dx = g->x[0] - potential->x[0];
+  const float dy = g->x[1] - potential->x[1];
+  const float dz = g->x[2] - potential->x[2];
 
   return -0.5f * potential->vrot * potential->vrot *
          logf(dx * dx + dy * dy + dz * dz + potential->epsilon2);
@@ -166,15 +166,12 @@ static INLINE void potential_init_backend(
     const struct unit_system* us, const struct space* s,
     struct external_potential* potential) {
 
-  potential->x =
-      s->dim[0] / 2. +
-      parser_get_param_double(parameter_file, "IsothermalPotential:position_x");
-  potential->y =
-      s->dim[1] / 2. +
-      parser_get_param_double(parameter_file, "IsothermalPotential:position_y");
-  potential->z =
-      s->dim[2] / 2. +
-      parser_get_param_double(parameter_file, "IsothermalPotential:position_z");
+  parser_get_param_double_array(parameter_file, "IsothermalPotential:position",
+                                3, potential->x);
+  potential->x[0] += s->dim[0] / 2.;
+  potential->x[1] += s->dim[1] / 2.;
+  potential->x[2] += s->dim[2] / 2.;
+
   potential->vrot =
       parser_get_param_double(parameter_file, "IsothermalPotential:vrot");
   potential->timestep_mult = parser_get_param_float(
@@ -198,7 +195,7 @@ static INLINE void potential_print_backend(
       "External potential is 'Isothermal' with properties are (x,y,z) = (%e, "
       "%e, %e), vrot = %e "
       "timestep multiplier = %e, epsilon = %e",
-      potential->x, potential->y, potential->z, potential->vrot,
+      potential->x[0], potential->x[1], potential->x[2], potential->vrot,
       potential->timestep_mult, sqrtf(potential->epsilon2));
 }
 
diff --git a/src/potential/point_mass/potential.h b/src/potential/point_mass/potential.h
index db875842d51f6c0a28dd1308fd7dd1728e746ce4..f9d56a1ff165f2331c91ea828b5ffe0e0db76c2f 100644
--- a/src/potential/point_mass/potential.h
+++ b/src/potential/point_mass/potential.h
@@ -41,7 +41,7 @@
 struct external_potential {
 
   /*! Position of the point mass */
-  double x, y, z;
+  double x[3];
 
   /*! Mass */
   double mass;
@@ -66,15 +66,15 @@ __attribute__((always_inline)) INLINE static float external_gravity_timestep(
     const struct gpart* restrict g) {
 
   const float G_newton = phys_const->const_newton_G;
-  const float dx = g->x[0] - potential->x;
-  const float dy = g->x[1] - potential->y;
-  const float dz = g->x[2] - potential->z;
+  const float dx = g->x[0] - potential->x[0];
+  const float dy = g->x[1] - potential->x[1];
+  const float dz = g->x[2] - potential->x[2];
   const float rinv = 1.f / sqrtf(dx * dx + dy * dy + dz * dz);
   const float rinv2 = rinv * rinv;
   const float rinv3 = rinv2 * rinv;
-  const float drdv = (g->x[0] - potential->x) * (g->v_full[0]) +
-                     (g->x[1] - potential->y) * (g->v_full[1]) +
-                     (g->x[2] - potential->z) * (g->v_full[2]);
+  const float drdv = (g->x[0] - potential->x[0]) * (g->v_full[0]) +
+                     (g->x[1] - potential->x[1]) * (g->v_full[1]) +
+                     (g->x[2] - potential->x[2]) * (g->v_full[2]);
   const float dota_x = G_newton * potential->mass * rinv3 *
                        (-g->v_full[0] + 3.f * rinv2 * drdv * dx);
   const float dota_y = G_newton * potential->mass * rinv3 *
@@ -109,9 +109,9 @@ __attribute__((always_inline)) INLINE static void external_gravity_acceleration(
     double time, const struct external_potential* restrict potential,
     const struct phys_const* restrict phys_const, struct gpart* restrict g) {
 
-  const float dx = g->x[0] - potential->x;
-  const float dy = g->x[1] - potential->y;
-  const float dz = g->x[2] - potential->z;
+  const float dx = g->x[0] - potential->x[0];
+  const float dy = g->x[1] - potential->x[1];
+  const float dz = g->x[2] - potential->x[2];
   const float rinv = 1.f / sqrtf(dx * dx + dy * dy + dz * dz);
   const float rinv3 = rinv * rinv * rinv;
 
@@ -134,9 +134,9 @@ external_gravity_get_potential_energy(
     double time, const struct external_potential* potential,
     const struct phys_const* const phys_const, const struct gpart* g) {
 
-  const float dx = g->x[0] - potential->x;
-  const float dy = g->x[1] - potential->y;
-  const float dz = g->x[2] - potential->z;
+  const float dx = g->x[0] - potential->x[0];
+  const float dy = g->x[1] - potential->x[1];
+  const float dz = g->x[2] - potential->x[2];
   const float rinv = 1. / sqrtf(dx * dx + dy * dy + dz * dz);
   return -phys_const->const_newton_G * potential->mass * rinv;
 }
@@ -156,12 +156,8 @@ static INLINE void potential_init_backend(
     const struct unit_system* us, const struct space* s,
     struct external_potential* potential) {
 
-  potential->x =
-      parser_get_param_double(parameter_file, "PointMassPotential:position_x");
-  potential->y =
-      parser_get_param_double(parameter_file, "PointMassPotential:position_y");
-  potential->z =
-      parser_get_param_double(parameter_file, "PointMassPotential:position_z");
+  parser_get_param_double_array(parameter_file, "PointMassPotential:position",
+                                3, potential->x);
   potential->mass =
       parser_get_param_double(parameter_file, "PointMassPotential:mass");
   potential->timestep_mult = parser_get_param_float(
@@ -179,7 +175,7 @@ static INLINE void potential_print_backend(
   message(
       "External potential is 'Point mass' with properties (x,y,z) = (%e, %e, "
       "%e), M = %e timestep multiplier = %e.",
-      potential->x, potential->y, potential->z, potential->mass,
+      potential->x[0], potential->x[1], potential->x[2], potential->mass,
       potential->timestep_mult);
 }
 
diff --git a/src/potential/point_mass_ring/potential.h b/src/potential/point_mass_ring/potential.h
index 551efe32521a5c5ee8068ba409dbb81547103e8f..b384c3df6771c9f3342f5aed6514e75c77bf2c17 100644
--- a/src/potential/point_mass_ring/potential.h
+++ b/src/potential/point_mass_ring/potential.h
@@ -41,7 +41,7 @@
 struct external_potential {
 
   /*! Position of the point mass */
-  double x, y, z;
+  double x[3];
 
   /*! Mass */
   double mass;
@@ -68,16 +68,16 @@ __attribute__((always_inline)) INLINE static float external_gravity_timestep(
     const struct gpart* restrict g) {
 
   const float G_newton = phys_const->const_newton_G;
-  const float dx = g->x[0] - potential->x;
-  const float dy = g->x[1] - potential->y;
-  const float dz = g->x[2] - potential->z;
+  const float dx = g->x[0] - potential->x[0];
+  const float dy = g->x[1] - potential->x[1];
+  const float dz = g->x[2] - potential->x[2];
   const float r = sqrtf(dx * dx + dy * dy + dz * dz);
   const float rinv = 1.f / r;
   const float rinv2 = rinv * rinv;
   const float rinv3 = rinv2 * rinv;
-  const float drdv = (g->x[0] - potential->x) * (g->v_full[0]) +
-                     (g->x[1] - potential->y) * (g->v_full[1]) +
-                     (g->x[2] - potential->z) * (g->v_full[2]);
+  const float drdv = (g->x[0] - potential->x[0]) * (g->v_full[0]) +
+                     (g->x[1] - potential->x[1]) * (g->v_full[1]) +
+                     (g->x[2] - potential->x[2]) * (g->v_full[2]);
   float factor;
 
   if (r < 0.175) {
@@ -129,9 +129,9 @@ __attribute__((always_inline)) INLINE static void external_gravity_acceleration(
     double time, const struct external_potential* restrict potential,
     const struct phys_const* restrict phys_const, struct gpart* restrict g) {
 
-  const float dx = g->x[0] - potential->x;
-  const float dy = g->x[1] - potential->y;
-  const float dz = g->x[2] - potential->z;
+  const float dx = g->x[0] - potential->x[0];
+  const float dy = g->x[1] - potential->x[1];
+  const float dz = g->x[2] - potential->x[2];
   const float r = sqrtf(dx * dx + dy * dy + dz * dz);
   const float rinv = 1.f / r;
   const float rinv2 = rinv * rinv;
@@ -174,9 +174,9 @@ external_gravity_get_potential_energy(
     double time, const struct external_potential* potential,
     const struct phys_const* const phys_const, const struct gpart* g) {
 
-  const float dx = g->x[0] - potential->x;
-  const float dy = g->x[1] - potential->y;
-  const float dz = g->x[2] - potential->z;
+  const float dx = g->x[0] - potential->x[0];
+  const float dy = g->x[1] - potential->x[1];
+  const float dz = g->x[2] - potential->x[2];
   const float rinv = 1. / sqrtf(dx * dx + dy * dy + dz * dz);
   return -phys_const->const_newton_G * potential->mass * rinv;
 }
@@ -196,12 +196,8 @@ static INLINE void potential_init_backend(
     const struct unit_system* us, const struct space* s,
     struct external_potential* potential) {
 
-  potential->x =
-      parser_get_param_double(parameter_file, "PointMassPotential:position_x");
-  potential->y =
-      parser_get_param_double(parameter_file, "PointMassPotential:position_y");
-  potential->z =
-      parser_get_param_double(parameter_file, "PointMassPotential:position_z");
+  parser_get_param_double_array(parameter_file, "PointMassPotential:position",
+                                3, potential->x);
   potential->mass =
       parser_get_param_double(parameter_file, "PointMassPotential:mass");
   potential->timestep_mult = parser_get_param_float(
@@ -219,7 +215,7 @@ static INLINE void potential_print_backend(
   message(
       "External potential is 'Point mass' with properties (x,y,z) = (%e, %e, "
       "%e), M = %e timestep multiplier = %e.",
-      potential->x, potential->y, potential->z, potential->mass,
+      potential->x[0], potential->x[1], potential->x[2], potential->mass,
       potential->timestep_mult);
 }
 
diff --git a/src/potential/point_mass_softened/potential.h b/src/potential/point_mass_softened/potential.h
index 80959ec923cbedcbea5fba3293c8ae4f94f65679..0e35e7bb9870c7954b47316a3cc30bb68cde5fc4 100644
--- a/src/potential/point_mass_softened/potential.h
+++ b/src/potential/point_mass_softened/potential.h
@@ -50,7 +50,7 @@
 struct external_potential {
 
   /*! Position of the point mass */
-  double x, y, z;
+  double x[3];
 
   /*! Mass */
   double mass;
@@ -77,9 +77,9 @@ __attribute__((always_inline)) INLINE static float external_gravity_timestep(
     const struct phys_const* restrict phys_const,
     const struct gpart* restrict g) {
 
-  const float dx = g->x[0] - potential->x;
-  const float dy = g->x[1] - potential->y;
-  const float dz = g->x[2] - potential->z;
+  const float dx = g->x[0] - potential->x[0];
+  const float dy = g->x[1] - potential->x[1];
+  const float dz = g->x[2] - potential->x[2];
 
   const float softening2 = potential->softening * potential->softening;
   const float r2 = dx * dx + dy * dy + dz * dz;
@@ -133,9 +133,9 @@ __attribute__((always_inline)) INLINE static void external_gravity_acceleration(
     double time, const struct external_potential* restrict potential,
     const struct phys_const* restrict phys_const, struct gpart* restrict g) {
 
-  const float dx = g->x[0] - potential->x;
-  const float dy = g->x[1] - potential->y;
-  const float dz = g->x[2] - potential->z;
+  const float dx = g->x[0] - potential->x[0];
+  const float dy = g->x[1] - potential->x[1];
+  const float dz = g->x[2] - potential->x[2];
   const float rinv = 1.f / sqrtf(dx * dx + dy * dy + dz * dz +
                                  potential->softening * potential->softening);
   const float rinv3 = rinv * rinv * rinv;
@@ -159,9 +159,9 @@ external_gravity_get_potential_energy(
     double time, const struct external_potential* potential,
     const struct phys_const* const phys_const, const struct gpart* g) {
 
-  const float dx = g->x[0] - potential->x;
-  const float dy = g->x[1] - potential->y;
-  const float dz = g->x[2] - potential->z;
+  const float dx = g->x[0] - potential->x[0];
+  const float dy = g->x[1] - potential->x[1];
+  const float dz = g->x[2] - potential->x[2];
   const float rinv = 1. / sqrtf(dx * dx + dy * dy + dz * dz +
                                 potential->softening * potential->softening);
 
@@ -183,12 +183,8 @@ static INLINE void potential_init_backend(
     const struct unit_system* us, const struct space* s,
     struct external_potential* potential) {
 
-  potential->x =
-      parser_get_param_double(parameter_file, "PointMassPotential:position_x");
-  potential->y =
-      parser_get_param_double(parameter_file, "PointMassPotential:position_y");
-  potential->z =
-      parser_get_param_double(parameter_file, "PointMassPotential:position_z");
+  parser_get_param_double_array(parameter_file, "PointMassPotential:position",
+                                3, potential->x);
   potential->mass =
       parser_get_param_double(parameter_file, "PointMassPotential:mass");
   potential->timestep_mult = parser_get_param_float(
@@ -208,7 +204,7 @@ static INLINE void potential_print_backend(
   message(
       "External potential is 'Point mass' with properties (x,y,z) = (%e, %e, "
       "%e), M = %e timestep multiplier = %e, softening = %e.",
-      potential->x, potential->y, potential->z, potential->mass,
+      potential->x[0], potential->x[1], potential->x[2], potential->mass,
       potential->timestep_mult, potential->softening);
 }
 
diff --git a/src/space.c b/src/space.c
index 5072318e690010973210a6a905c64525b7c4cce2..9fec6cb5e3f094716850f00171504145985713e9 100644
--- a/src/space.c
+++ b/src/space.c
@@ -2758,12 +2758,8 @@ void space_init(struct space *s, struct swift_params *params,
 
   /* Apply shift */
   double shift[3] = {0.0, 0.0, 0.0};
-  shift[0] =
-      parser_get_opt_param_double(params, "InitialConditions:shift_x", 0.0);
-  shift[1] =
-      parser_get_opt_param_double(params, "InitialConditions:shift_y", 0.0);
-  shift[2] =
-      parser_get_opt_param_double(params, "InitialConditions:shift_z", 0.0);
+  parser_get_opt_param_double_array(params, "InitialConditions:shift", 3,
+                                    shift);
   if ((shift[0] != 0. || shift[1] != 0. || shift[2] != 0.) && !dry_run) {
     message("Shifting particles by [%e %e %e]", shift[0], shift[1], shift[2]);
     for (size_t k = 0; k < Npart; k++) {