diff --git a/Documentation/ProjectFile/prj/time_loop/output/hdf/i_hdf.md b/Documentation/ProjectFile/prj/time_loop/output/hdf/i_hdf.md
new file mode 100644
index 0000000000000000000000000000000000000000..2043711a7a77bcb4ca95a5845138c5bfb9b6761d
--- /dev/null
+++ b/Documentation/ProjectFile/prj/time_loop/output/hdf/i_hdf.md
@@ -0,0 +1,2 @@
+Group of parameters when XDMF/HDF writer is used.
+Type of output in time_loop must be XDMF
diff --git a/Documentation/ProjectFile/prj/time_loop/output/hdf/t_number_of_files.md b/Documentation/ProjectFile/prj/time_loop/output/hdf/t_number_of_files.md
new file mode 100644
index 0000000000000000000000000000000000000000..662f4178c831d66745976513dc56e3a36f240376
--- /dev/null
+++ b/Documentation/ProjectFile/prj/time_loop/output/hdf/t_number_of_files.md
@@ -0,0 +1 @@
+\copydoc ProcessLib::Output::_n_files
diff --git a/MeshLib/IO/XDMF/HdfData.cpp b/MeshLib/IO/XDMF/HdfData.cpp
index 76b3ce2da9a83ebb5fec5a0b5361317ede20ae23..dc0b6272af6f1332a5375ed02ec7b2cfd4ec1eee 100644
--- a/MeshLib/IO/XDMF/HdfData.cpp
+++ b/MeshLib/IO/XDMF/HdfData.cpp
@@ -42,10 +42,12 @@ static hid_t meshPropertyType2HdfType(MeshPropertyDataType const ogs_data_type)
 
 HdfData::HdfData(void const* data_start, std::size_t const size_partitioned_dim,
                  std::size_t const size_tuple, std::string const& name,
-                 MeshPropertyDataType const mesh_property_data_type)
+                 MeshPropertyDataType const mesh_property_data_type,
+                 unsigned int const num_of_files)
     : data_start(data_start), name(name)
 {
-    auto const& partition_info = getPartitionInfo(size_partitioned_dim);
+    auto const& partition_info =
+        getPartitionInfo(size_partitioned_dim, num_of_files);
     auto const& offset_partitioned_dim = partition_info.local_offset;
     offsets = {offset_partitioned_dim, 0};
 
diff --git a/MeshLib/IO/XDMF/HdfData.h b/MeshLib/IO/XDMF/HdfData.h
index 50f6eb17a5be138b362bf03d8b5ddac378b0c0a5..bc41db110b2f0d10a0b5b0fcb9b409c3f2d72f82 100644
--- a/MeshLib/IO/XDMF/HdfData.h
+++ b/MeshLib/IO/XDMF/HdfData.h
@@ -25,7 +25,8 @@ struct HdfData final
 {
     HdfData(void const* data_start, std::size_t size_partitioned_dim,
             std::size_t size_tuple, std::string const& name,
-            MeshPropertyDataType mesh_property_data_type);
+            MeshPropertyDataType mesh_property_data_type,
+            unsigned int num_of_files);
     void const* data_start;
     std::vector<Hdf5DimType> data_space;
     std::vector<Hdf5DimType> offsets;
diff --git a/MeshLib/IO/XDMF/HdfWriter.cpp b/MeshLib/IO/XDMF/HdfWriter.cpp
index 4723120c8a5b8d871fe8d94644bd0531cc5894af..a4ba35cb9c8f6b7bb606396ac49799107c502464 100644
--- a/MeshLib/IO/XDMF/HdfWriter.cpp
+++ b/MeshLib/IO/XDMF/HdfWriter.cpp
@@ -213,9 +213,10 @@ HdfWriter::HdfWriter(std::vector<MeshHdfData> meshes,
                      unsigned long long const initial_step,
                      std::filesystem::path const& filepath,
                      bool const use_compression,
-                     bool const is_file_manager)
+                     bool const is_file_manager,
+                     unsigned int const num_of_files)
     : _hdf5_filepath(filepath),
-      _file(createFile(filepath)),
+      _file(createFile(filepath, num_of_files)),
       _meshes_group(
           H5Gcreate2(_file, "/meshes", H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)),
       _step_times{0},  // ToDo need to be initial time
diff --git a/MeshLib/IO/XDMF/HdfWriter.h b/MeshLib/IO/XDMF/HdfWriter.h
index 97c8707850a552134f3efd4092502089933e88ee..a58c1f4bd4dca35ac1ee4f5033112d231418ddc2 100644
--- a/MeshLib/IO/XDMF/HdfWriter.h
+++ b/MeshLib/IO/XDMF/HdfWriter.h
@@ -43,13 +43,15 @@ public:
      * @param filepath absolute or relative filepath to the hdf5 file
      * @param use_compression if true gzip compression is enabled
      * @param is_file_manager True if process (in parallel execution) is
+     * @param num_of_files Number of outputfiles
      * File_Manager
      */
     HdfWriter(std::vector<MeshHdfData> meshes,
               unsigned long long initial_step,
               std::filesystem::path const& filepath,
               bool use_compression,
-              bool is_file_manager);
+              bool is_file_manager,
+              unsigned int num_of_files);
     /**
      * \brief Writes attributes. The data
      * itself is hold by a structure outside of this class. The writer assumes
diff --git a/MeshLib/IO/XDMF/XdmfData.cpp b/MeshLib/IO/XDMF/XdmfData.cpp
index 432895df556c906536ff1bd612fdfb0065b8eb04..f2e2b3c9eef39c19e8153cb39e652e434b02f5e2 100644
--- a/MeshLib/IO/XDMF/XdmfData.cpp
+++ b/MeshLib/IO/XDMF/XdmfData.cpp
@@ -23,7 +23,8 @@ XdmfData::XdmfData(std::size_t const size_partitioned_dim,
                    MeshPropertyDataType const mesh_property_data_type,
                    std::string const& name,
                    std::optional<MeshLib::MeshItemType> const attribute_center,
-                   unsigned int const index)
+                   unsigned int const index,
+                   unsigned int const num_of_files)
     : starts(
           [&size_tuple]()
           {
@@ -53,7 +54,7 @@ XdmfData::XdmfData(std::size_t const size_partitioned_dim,
       attribute_center(attribute_center),
       index(index)
 {
-    auto partition_info = getPartitionInfo(size_partitioned_dim);
+    auto partition_info = getPartitionInfo(size_partitioned_dim, num_of_files);
     // TODO (tm) XdmfLib does not support 64 bit data types so far
     assert(partition_info.local_length <
            std::numeric_limits<unsigned int>::max());
diff --git a/MeshLib/IO/XDMF/XdmfData.h b/MeshLib/IO/XDMF/XdmfData.h
index 5d256f248120b5ec1c44d2f5839b0e824126164e..d3a96f11d261c8294273b24ad06735195e3f03bc 100644
--- a/MeshLib/IO/XDMF/XdmfData.h
+++ b/MeshLib/IO/XDMF/XdmfData.h
@@ -46,13 +46,15 @@ struct XdmfData final
      * @param index The position of the DataItem parents in a grid
      * (representing a single step). Convention is: 1=Time, 2=
      * Geometry, 3=Topology, 4>=Attribute
+     * @param num_of_files If greater than 1 it groups the data of each process.
+     * The num_of_files specifies the number of groups
      *
      */
     XdmfData(std::size_t size_partitioned_dim, std::size_t size_tuple,
              MeshPropertyDataType mesh_property_data_type,
              std::string const& name,
              std::optional<MeshLib::MeshItemType> attribute_center,
-             unsigned int const index);
+             unsigned int const index, unsigned int num_of_files);
     // a hyperslab is defined by starts and strides see
     // https://www.xdmf.org/index.php/XDMF_Model_and_Format#HyperSlab
     std::vector<XdmfDimType> starts;
diff --git a/MeshLib/IO/XDMF/XdmfHdfWriter.cpp b/MeshLib/IO/XDMF/XdmfHdfWriter.cpp
index 3d3e1381d1606a2dd14105044aea8cc0cfb25910..265a639422550ce3f4938d67c8dd7ba2894e9d12 100644
--- a/MeshLib/IO/XDMF/XdmfHdfWriter.cpp
+++ b/MeshLib/IO/XDMF/XdmfHdfWriter.cpp
@@ -39,7 +39,7 @@ XdmfHdfWriter::XdmfHdfWriter(
     std::filesystem::path const& filepath, unsigned long long const time_step,
     double const initial_time,
     std::set<std::string> const& variable_output_names,
-    bool const use_compression)
+    bool const use_compression, unsigned int const num_of_files)
 {
     // ogs meshes to vector of Xdmf/HDF meshes (we keep Xdmf and HDF together
     // because XDMF depends on HDF) to meta
@@ -79,12 +79,12 @@ XdmfHdfWriter::XdmfHdfWriter(
     // Transform the data to be written into a format conforming with the rules
     // of xdmf topology and geometry
     auto const transform_ogs_mesh_data_to_xdmf_conforming_data =
-        [](auto const& mesh)
+        [&num_of_files](auto const& mesh)
     {
         auto flattened_geometry_values = transformToXDMFGeometry(mesh);
         // actually this line is only needed to calculate the offset
-        XdmfHdfData const& geometry =
-            transformGeometry(mesh, flattened_geometry_values.data());
+        XdmfHdfData const& geometry = transformGeometry(
+            mesh, flattened_geometry_values.data(), num_of_files);
         auto const flattened_topology_values =
             transformToXDMFTopology(mesh, geometry.hdf.offsets[0]);
         return std::make_unique<TransformedMeshData>(
@@ -94,17 +94,19 @@ XdmfHdfWriter::XdmfHdfWriter(
 
     // create metadata for transformed data and original ogs mesh data
     auto const transform_to_meta_data =
-        [&transform_ogs_mesh_data_to_xdmf_conforming_data](auto const& mesh)
+        [&transform_ogs_mesh_data_to_xdmf_conforming_data,
+         &num_of_files](auto const& mesh)
     {
         // important: transformed data must survive and be unique, raw pointer
         // to its memory!
         std::unique_ptr<TransformedMeshData> xdmf_conforming_data =
             transform_ogs_mesh_data_to_xdmf_conforming_data(mesh);
         auto const geometry = transformGeometry(
-            mesh, xdmf_conforming_data->flattened_geometry_values.data());
-        auto const topology =
-            transformTopology(xdmf_conforming_data->flattened_topology_values);
-        auto const attributes = transformAttributes(mesh);
+            mesh, xdmf_conforming_data->flattened_geometry_values.data(),
+            num_of_files);
+        auto const topology = transformTopology(
+            xdmf_conforming_data->flattened_topology_values, num_of_files);
+        auto const attributes = transformAttributes(mesh, num_of_files);
         return XdmfHdfMesh{std::move(geometry), std::move(topology),
                            std::move(attributes), mesh.get().getName(),
                            std::move(xdmf_conforming_data)};
@@ -157,7 +159,7 @@ XdmfHdfWriter::XdmfHdfWriter(
     auto const is_file_manager = isFileManager();
     _hdf_writer = std::make_unique<HdfWriter>(std::move(hdf_meshes), time_step,
                                               hdf_filepath, use_compression,
-                                              is_file_manager);
+                                              is_file_manager, num_of_files);
 
     // --------------- XDMF ---------------------
     // The light data is only written by just one process
@@ -174,7 +176,8 @@ XdmfHdfWriter::XdmfHdfWriter(
     {
         std::string const xdmf_name = metamesh.name;
         std::filesystem::path const xdmf_filepath =
-            filepath.parent_path() / (xdmf_name + ".xdmf");
+            filepath.parent_path() /
+            (filepath.stem().string() + "_" + xdmf_name + ".xdmf");
 
         std::vector<XdmfData> xdmf_attributes;
         std::transform(metamesh.attributes.begin(), metamesh.attributes.end(),
diff --git a/MeshLib/IO/XDMF/XdmfHdfWriter.h b/MeshLib/IO/XDMF/XdmfHdfWriter.h
index 69d8b8e233fbaff0327339a7ef4bf85300876506..bdb3b884ec3c1bd733dc6421422095ce09712c6c 100644
--- a/MeshLib/IO/XDMF/XdmfHdfWriter.h
+++ b/MeshLib/IO/XDMF/XdmfHdfWriter.h
@@ -36,12 +36,13 @@ public:
      * that change over time
      * @param use_compression if true, zlib compression in HDFWriter component
      * is used
+     * @param num_of_files number of hdf5 output files
      */
     XdmfHdfWriter(
         std::vector<std::reference_wrapper<const MeshLib::Mesh>> meshes,
         std::filesystem::path const& filepath, unsigned long long time_step,
         double initial_time, std::set<std::string> const& variable_output_names,
-        bool use_compression);
+        bool use_compression, unsigned int num_of_files);
 
     /**
      * \brief Adds data for either lazy (xdmf) or eager (hdf) writing algorithm
diff --git a/MeshLib/IO/XDMF/fileIO.h b/MeshLib/IO/XDMF/fileIO.h
index 1fe20cd1e448ec02edd7133f43fa6c67445a95a9..1e03fafa175c1b474b831faed32ca5ad24914925 100644
--- a/MeshLib/IO/XDMF/fileIO.h
+++ b/MeshLib/IO/XDMF/fileIO.h
@@ -18,7 +18,9 @@
 
 namespace MeshLib::IO
 {
-int64_t createFile(std::filesystem::path const& filepath);
-int64_t openHDF5File(std::filesystem::path const& filepath);
+int64_t createFile(std::filesystem::path const& filepath,
+                   unsigned int num_of_files);
+int64_t openHDF5File(std::filesystem::path const& filepath,
+                     unsigned int num_of_files);
 int64_t createHDF5TransferPolicy();
 }  // namespace MeshLib::IO
diff --git a/MeshLib/IO/XDMF/mpi/fileIO.cpp b/MeshLib/IO/XDMF/mpi/fileIO.cpp
index 6b3215f61c10c9102097272b97a4f808a7fac051..6e1343056508108386a98ec75991c00d82cd470b 100644
--- a/MeshLib/IO/XDMF/mpi/fileIO.cpp
+++ b/MeshLib/IO/XDMF/mpi/fileIO.cpp
@@ -17,24 +17,47 @@
 #include <mpi.h>
 
 #include "BaseLib/Logging.h"
+#include "getCommunicator.h"
 
+using namespace std::string_literals;
 namespace MeshLib::IO
 {
-hid_t createFile(std::filesystem::path const& filepath)
+std::filesystem::path partitionFilename(
+    std::filesystem::path const& basic_filepath, int file_group)
 {
-    MPI_Comm comm = MPI_COMM_WORLD;
-    MPI_Info info = MPI_INFO_NULL;
+    std::string const filename = (file_group > 0)
+                                     ? basic_filepath.stem().string() + "_"s +
+                                           std::to_string(file_group) +
+                                           basic_filepath.extension().string()
+                                     : basic_filepath.filename().string();
+    std::filesystem::path const filepathwithextension =
+        basic_filepath.parent_path() / filename;
+    DBUG("HDF Filepath: {:s}.", filepathwithextension.string());
+    return filepathwithextension;
+};
+
+hid_t createFile(std::filesystem::path const& filepath,
+                 unsigned int const num_of_files)
+{
+    auto const communicator = getCommunicator(num_of_files);
+    MPI_Comm const comm = communicator.mpi_communicator;
+    MPI_Info const info = MPI_INFO_NULL;
     hid_t const plist_id = H5Pcreate(H5P_FILE_ACCESS);
+
     H5Pset_fapl_mpio(plist_id, comm, info);
-    hid_t file = H5Fcreate(filepath.string().c_str(), H5F_ACC_TRUNC,
+    std::filesystem::path const partition_filename =
+        partitionFilename(filepath, communicator.color);
+    hid_t file = H5Fcreate(partition_filename.string().c_str(), H5F_ACC_TRUNC,
                            H5P_DEFAULT, plist_id);
     H5Pclose(plist_id);
+
     return file;
 }
 
-hid_t openHDF5File(std::filesystem::path const& filepath)
+hid_t openHDF5File(std::filesystem::path const& filepath,
+                   unsigned int const num_of_files)
 {
-    MPI_Comm comm = MPI_COMM_WORLD;
+    MPI_Comm const comm = getCommunicator(num_of_files).mpi_communicator;
     MPI_Info info = MPI_INFO_NULL;
     hid_t const plist_id = H5Pcreate(H5P_FILE_ACCESS);
     H5Pset_fapl_mpio(plist_id, comm, info);
diff --git a/MeshLib/IO/XDMF/mpi/getCommunicator.cpp b/MeshLib/IO/XDMF/mpi/getCommunicator.cpp
new file mode 100644
index 0000000000000000000000000000000000000000..c6ac67395e772e63d9c29308e7be8336a80f759c
--- /dev/null
+++ b/MeshLib/IO/XDMF/mpi/getCommunicator.cpp
@@ -0,0 +1,55 @@
+/**
+ * \file
+ * \author Tobias Meisel
+ * \date   2020-12-08
+ * \brief  Function specific to execution with MPI, never include directly!!
+ * \copyright
+ * Copyright (c) 2012-2021, OpenGeoSys Community (http://www.opengeosys.org)
+ *            Distributed under a Modified BSD License.
+ *              See accompanying file LICENSE.txt or
+ *              http://www.opengeosys.org/project/license
+ *
+ */
+
+#include "getCommunicator.h"
+
+#include <hdf5.h>
+#include <math.h>
+#include <mpi.h>
+
+#include <cassert>
+
+#include "BaseLib/Logging.h"
+#include "MeshLib/IO/XDMF/fileIO.h"
+
+using namespace std::string_literals;
+
+namespace MeshLib::IO
+{
+int getGroupIndex(int const input_index, int const input_size,
+                  int const new_group_size)
+{
+    // A grouping algorithm that determines the number of groups and return the
+    // group idx of the specified input_index
+    assert(input_size >= new_group_size);
+    int const minimal_output_group_size =
+        std::lround(input_size / new_group_size);
+    int const maximal_output_group_size = (input_size % new_group_size)
+                                              ? minimal_output_group_size + 1
+                                              : minimal_output_group_size;
+    return std::lround(input_index / maximal_output_group_size);
+};
+
+FileCommunicator getCommunicator(unsigned int const num_of_files)
+{
+    int num_procs;
+    MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
+    int rank_id;
+    MPI_Comm_rank(MPI_COMM_WORLD, &rank_id);
+    int file_group_id = getGroupIndex(rank_id, num_procs, num_of_files);
+    MPI_Comm new_communicator;
+    MPI_Comm_split(MPI_COMM_WORLD, file_group_id, rank_id, &new_communicator);
+    return FileCommunicator{std::move(new_communicator),
+                            std::move(file_group_id), ""};
+}
+}  // namespace MeshLib::IO
\ No newline at end of file
diff --git a/MeshLib/IO/XDMF/mpi/getCommunicator.h b/MeshLib/IO/XDMF/mpi/getCommunicator.h
new file mode 100644
index 0000000000000000000000000000000000000000..416b59ec2a62e0f0c40703bca0e042ea8c67375b
--- /dev/null
+++ b/MeshLib/IO/XDMF/mpi/getCommunicator.h
@@ -0,0 +1,16 @@
+#pragma once
+
+#include <mpi.h>
+
+#include <filesystem>
+
+namespace MeshLib::IO
+{
+struct FileCommunicator final
+{
+    MPI_Comm mpi_communicator;
+    int color;
+    std::filesystem::path output_filename;
+};
+FileCommunicator getCommunicator(unsigned int num_of_files);
+}  // namespace MeshLib::IO
\ No newline at end of file
diff --git a/MeshLib/IO/XDMF/mpi/partition.cpp b/MeshLib/IO/XDMF/mpi/partition.cpp
index a398f6a5f383d3d3fb75b65723bc2fa148fe535e..67e0e5791faa1d01260ce9e9711cd704261f5664 100644
--- a/MeshLib/IO/XDMF/mpi/partition.cpp
+++ b/MeshLib/IO/XDMF/mpi/partition.cpp
@@ -15,11 +15,11 @@
 
 #include <mpi.h>
 
-#include <deque>
 #include <numeric>
 
 #include "BaseLib/Logging.h"
 #include "MeshLib/IO/XDMF/fileIO.h"
+#include "getCommunicator.h"
 
 namespace MeshLib::IO
 {
@@ -30,9 +30,10 @@ bool isFileManager()
     return mpi_rank == 0;
 }
 
-PartitionInfo getPartitionInfo(std::size_t const size)
+PartitionInfo getPartitionInfo(std::size_t const size,
+                               unsigned int const num_of_files)
 {
-    MPI_Comm const mpi_comm = MPI_COMM_WORLD;
+    MPI_Comm const mpi_comm = getCommunicator(num_of_files).mpi_communicator;
     int mpi_size;
     int mpi_rank;
     MPI_Comm_size(mpi_comm, &mpi_size);
diff --git a/MeshLib/IO/XDMF/partition.h b/MeshLib/IO/XDMF/partition.h
index a343c590fe9d969397fe3790783a5b4fae4f6bae..afb2cb8134d91eda895581513af5ef5988c0a7de 100644
--- a/MeshLib/IO/XDMF/partition.h
+++ b/MeshLib/IO/XDMF/partition.h
@@ -24,6 +24,7 @@ struct PartitionInfo
     std::size_t global_length;
 };
 
-PartitionInfo getPartitionInfo(std::size_t const size);
+PartitionInfo getPartitionInfo(std::size_t const size,
+                               unsigned int const num_of_files);
 bool isFileManager();
 }  // namespace MeshLib::IO
diff --git a/MeshLib/IO/XDMF/posix/fileIO.cpp b/MeshLib/IO/XDMF/posix/fileIO.cpp
index 2f7d081a8417d82729a41a5e0bb161e593a1fd6c..a9c3ad6b6fdd083bade4be547e096007e463f260 100644
--- a/MeshLib/IO/XDMF/posix/fileIO.cpp
+++ b/MeshLib/IO/XDMF/posix/fileIO.cpp
@@ -16,7 +16,7 @@
 #include <hdf5.h>
 namespace MeshLib::IO
 {
-int64_t createFile(std::filesystem::path const& filepath)
+int64_t createFile(std::filesystem::path const& filepath, unsigned int)
 {
     return H5Fcreate(filepath.string().c_str(), H5F_ACC_TRUNC, H5P_DEFAULT,
                      H5P_DEFAULT);
diff --git a/MeshLib/IO/XDMF/posix/partition.cpp b/MeshLib/IO/XDMF/posix/partition.cpp
index 2fe296585719e16226f4358d4837835e82baff92..15796288e1e1cd06bad2ece80a36b9713ecce6d8 100644
--- a/MeshLib/IO/XDMF/posix/partition.cpp
+++ b/MeshLib/IO/XDMF/posix/partition.cpp
@@ -17,7 +17,7 @@ bool isFileManager()
     return true;
 }
 
-PartitionInfo getPartitionInfo(std::size_t const size)
+PartitionInfo getPartitionInfo(std::size_t const size, unsigned int)
 {
         // local_offset, local_length, longest_local_length, global_number_process
     return {0, size, size, size};
diff --git a/MeshLib/IO/XDMF/transformData.cpp b/MeshLib/IO/XDMF/transformData.cpp
index cd07b709a27ac7058db884c36d67e3d544e1ccf5..c2bb18dfee53cbd48a6d32567d71c5164207ac89 100644
--- a/MeshLib/IO/XDMF/transformData.cpp
+++ b/MeshLib/IO/XDMF/transformData.cpp
@@ -70,7 +70,8 @@ constexpr auto cellTypeOGS2XDMF(MeshLib::CellType const& cell_type)
 }
 
 std::optional<XdmfHdfData> transformAttribute(
-    std::pair<std::string, PropertyVectorBase*> const& property_pair)
+    std::pair<std::string, PropertyVectorBase*> const& property_pair,
+    unsigned int const num_of_files)
 {
     // 3 data that will be captured and written by lambda f below
     MeshPropertyDataType data_type = MeshPropertyDataType::unknown;
@@ -190,16 +191,18 @@ std::optional<XdmfHdfData> transformAttribute(
 
     std::string const& name = property_base->getPropertyName();
 
-    HdfData hdf = {data_ptr, num_of_tuples, ui_global_components, name,
-                   data_type};
+    HdfData hdf = {data_ptr, num_of_tuples, ui_global_components,
+                   name,     data_type,     num_of_files};
 
     XdmfData xdmf = {num_of_tuples, ui_global_components, data_type,
-                     name,          mesh_item_type,       0};
+                     name,          mesh_item_type,       0,
+                     num_of_files};
 
     return XdmfHdfData{std::move(hdf), std::move(xdmf)};
 }
 
-std::vector<XdmfHdfData> transformAttributes(MeshLib::Mesh const& mesh)
+std::vector<XdmfHdfData> transformAttributes(MeshLib::Mesh const& mesh,
+                                             unsigned int const num_of_files)
 {
     MeshLib::Properties const& properties = mesh.getProperties();
 
@@ -213,8 +216,8 @@ std::vector<XdmfHdfData> transformAttributes(MeshLib::Mesh const& mesh)
             continue;
         }
 
-        if (auto const attribute =
-                transformAttribute(std::pair(name, property_base)))
+        if (auto const attribute = transformAttribute(
+                std::pair(name, property_base), num_of_files))
         {
             attributes.push_back(attribute.value());
         }
@@ -242,7 +245,9 @@ std::vector<double> transformToXDMFGeometry(MeshLib::Mesh const& mesh)
     return values;
 }
 
-XdmfHdfData transformGeometry(MeshLib::Mesh const& mesh, double const* data_ptr)
+XdmfHdfData transformGeometry(MeshLib::Mesh const& mesh,
+                              double const* data_ptr,
+                              unsigned int const num_of_files)
 {
     std::string const name = "geometry";
     std::vector<MeshLib::Node*> const& nodes = mesh.getNodes();
@@ -250,11 +255,16 @@ XdmfHdfData transformGeometry(MeshLib::Mesh const& mesh, double const* data_ptr)
     int const point_size = 3;
     auto const& partition_dim = nodes.size();
 
-    HdfData const hdf = {data_ptr, partition_dim, point_size, name,
-                         MeshPropertyDataType::float64};
+    HdfData const hdf = {data_ptr,
+                         partition_dim,
+                         point_size,
+                         name,
+                         MeshPropertyDataType::float64,
+                         num_of_files};
     XdmfData const xdmf = {
         partition_dim, point_size,   MeshPropertyDataType::float64,
-        name,          std::nullopt, 2};
+        name,          std::nullopt, 2,
+        num_of_files};
 
     return XdmfHdfData{std::move(hdf), std::move(xdmf)};
 }
@@ -286,13 +296,16 @@ std::vector<int> transformToXDMFTopology(MeshLib::Mesh const& mesh,
     return values;
 }
 
-XdmfHdfData transformTopology(std::vector<int> const& values)
+XdmfHdfData transformTopology(std::vector<int> const& values,
+                              unsigned int const num_of_files)
 {
     std::string const name = "topology";
-    HdfData const hdf = {values.data(), values.size(), 1, name,
-                         MeshPropertyDataType::int32};
-    XdmfData const xdmf = {values.size(), 1, MeshPropertyDataType::int32, name,
-                           std::nullopt,  3};
+    HdfData const hdf = {
+        values.data(), values.size(), 1, name, MeshPropertyDataType::int32,
+        num_of_files};
+    XdmfData const xdmf = {
+        values.size(), 1, MeshPropertyDataType::int32, name, std::nullopt, 3,
+        num_of_files};
 
     return XdmfHdfData{std::move(hdf), std::move(xdmf)};
 }
diff --git a/MeshLib/IO/XDMF/transformData.h b/MeshLib/IO/XDMF/transformData.h
index 16284421061aae9e1222594528f55b1344a47e0f..8c0e789cd5f839cde13b59c20dee588deb2caa1a 100644
--- a/MeshLib/IO/XDMF/transformData.h
+++ b/MeshLib/IO/XDMF/transformData.h
@@ -26,23 +26,31 @@ namespace MeshLib::IO
 /**
  * \brief Create meta data for attributes used for hdf5 and xdmf
  * @param mesh OGS mesh can be mesh or partitionedMesh
+ * @param num_of_files If greater than 1 it groups the data of each process. The
+ * num_of_files specifies the number of groups
  * @return vector of meta data
  */
-std::vector<XdmfHdfData> transformAttributes(MeshLib::Mesh const& mesh);
+std::vector<XdmfHdfData> transformAttributes(MeshLib::Mesh const& mesh,
+                                             unsigned int num_of_files);
 /**
  * \brief Create meta data for geometry used for hdf5 and xdmf
  * @param mesh OGS mesh can be mesh or partitionedMesh
  * @param data_ptr Memory location of geometry values.
+ * @param num_of_files If greater than 1 it groups the data of each process. The
+ * num_of_files specifies the number of groups
  * @return Geometry meta data
  */
-XdmfHdfData transformGeometry(MeshLib::Mesh const& mesh,
-                              double const* data_ptr);
+XdmfHdfData transformGeometry(MeshLib::Mesh const& mesh, double const* data_ptr,
+                              unsigned int num_of_files);
 /**
  * \brief  Create meta data for topology used for HDF5 and XDMF
  * @param values actual topology values to get size and memory location
+ * @param num_of_files If greater than 1 it groups the data of each process. The
+ * num_of_files specifies the number of groups
  * @return Topology meta data
  */
-XdmfHdfData transformTopology(std::vector<int> const& values);
+XdmfHdfData transformTopology(std::vector<int> const& values,
+                              unsigned int num_of_files);
 /**
  * \brief Copies all node points into a new vector. Contiguous data used for
  * writing. Conform with XDMF standard in
diff --git a/MeshLib/IO/writeMeshToFile.cpp b/MeshLib/IO/writeMeshToFile.cpp
index 36b0ef47815aee160fc0458641ff204eecfb9fa9..bf170fca3263d1fb1ffca50647422d9548c24135 100644
--- a/MeshLib/IO/writeMeshToFile.cpp
+++ b/MeshLib/IO/writeMeshToFile.cpp
@@ -51,7 +51,7 @@ int writeMeshToFile(const MeshLib::Mesh& mesh,
         const std::reference_wrapper<const MeshLib::Mesh> mr = mesh;
         meshes.push_back(mr);
         MeshLib::IO::XdmfHdfWriter(std::move(meshes), file_path, 0, 0.0,
-                                   variable_output_names, true);
+                                   variable_output_names, true, 1);
         return 0;
     }
     ERR("writeMeshToFile(): Unknown file extension '{:s}'. Can not write file "
diff --git a/ProcessLib/LiquidFlow/Tests.cmake b/ProcessLib/LiquidFlow/Tests.cmake
index 401fc4ae88a0f54188f46c0f8fde178e99fd8a45..5f06a0792dc390902be5d89b58966f72eb636c9f 100644
--- a/ProcessLib/LiquidFlow/Tests.cmake
+++ b/ProcessLib/LiquidFlow/Tests.cmake
@@ -483,16 +483,16 @@ AddTest(
     TESTER xdmfdiff
     REQUIREMENTS NOT OGS_USE_MPI
     DIFF_DATA
-    square_5x5_tris_32.xdmf square_5x5_tris_32.xdmf pressure pressure 1e-7 1e-13
-    square_5x5_tris_32.xdmf square_5x5_tris_32.xdmf HydraulicFlow HydraulicFlow 1e-7 1e-13
-    square_5x5_tris_32.xdmf square_5x5_tris_32.xdmf MaterialIDs MaterialIDs 1e-7 1e-13
-    square_5x5_tris_32.xdmf square_5x5_tris_32.xdmf v v 1e-7 1e-13
-    square_5x5_tris_32_right_boundary.xdmf square_5x5_tris_32_right_boundary.xdmf pressure pressure 1e-7 1e-13
-    square_5x5_tris_32_right_boundary.xdmf square_5x5_tris_32_right_boundary.xdmf bulk_element_ids bulk_element_ids 1e-7 1e-13
-    square_5x5_tris_32_right_boundary.xdmf square_5x5_tris_32_right_boundary.xdmf bulk_node_ids bulk_node_ids 1e-7 1e-13
-    square_5x5_tris_32_left_boundary.xdmf square_5x5_tris_32_left_boundary.xdmf pressure pressure 1e-7 1e-13
-    square_5x5_tris_32_left_boundary.xdmf square_5x5_tris_32_left_boundary.xdmf bulk_element_ids bulk_element_ids 1e-7 1e-13
-    square_5x5_tris_32_left_boundary.xdmf square_5x5_tris_32_left_boundary.xdmf bulk_node_ids bulk_node_ids 1e-7 1e-13
+    square_5x5_tris_32.xdmf square_5x5_tris_32_square_5x5_tris_32.xdmf pressure pressure 1e-7 1e-13
+    square_5x5_tris_32.xdmf square_5x5_tris_32_square_5x5_tris_32.xdmf HydraulicFlow HydraulicFlow 1e-7 1e-13
+    square_5x5_tris_32.xdmf square_5x5_tris_32_square_5x5_tris_32.xdmf MaterialIDs MaterialIDs 1e-7 1e-13
+    square_5x5_tris_32.xdmf square_5x5_tris_32_square_5x5_tris_32.xdmf v v 1e-7 1e-13
+    square_5x5_tris_32_right_boundary.xdmf square_5x5_tris_32_square_5x5_tris_32_right_boundary.xdmf pressure pressure 1e-7 1e-13
+    square_5x5_tris_32_right_boundary.xdmf square_5x5_tris_32_square_5x5_tris_32_right_boundary.xdmf bulk_element_ids bulk_element_ids 1e-7 1e-13
+    square_5x5_tris_32_right_boundary.xdmf square_5x5_tris_32_square_5x5_tris_32_right_boundary.xdmf bulk_node_ids bulk_node_ids 1e-7 1e-13
+    square_5x5_tris_32_left_boundary.xdmf square_5x5_tris_32_square_5x5_tris_32_left_boundary.xdmf pressure pressure 1e-7 1e-13
+    square_5x5_tris_32_left_boundary.xdmf square_5x5_tris_32_square_5x5_tris_32_left_boundary.xdmf bulk_element_ids bulk_element_ids 1e-7 1e-13
+    square_5x5_tris_32_left_boundary.xdmf square_5x5_tris_32_square_5x5_tris_32_left_boundary.xdmf bulk_node_ids bulk_node_ids 1e-7 1e-13
 )
 
 AddTest(
@@ -504,16 +504,16 @@ AddTest(
     TESTER xdmfdiff
     REQUIREMENTS NOT OGS_USE_MPI
     DIFF_DATA
-    square_5x5_tris_32.xdmf square_5x5_tris_32.xdmf pressure pressure 1e-7 1e-13
-    square_5x5_tris_32.xdmf square_5x5_tris_32.xdmf HydraulicFlow HydraulicFlow 1e-7 1e-13
-    square_5x5_tris_32.xdmf square_5x5_tris_32.xdmf MaterialIDs MaterialIDs 1e-7 1e-13
-    square_5x5_tris_32.xdmf square_5x5_tris_32.xdmf v v 1e-7 1e-13
-    square_5x5_tris_32_right_boundary.xdmf square_5x5_tris_32_right_boundary.xdmf pressure pressure 1e-7 1e-13
-    square_5x5_tris_32_right_boundary.xdmf square_5x5_tris_32_right_boundary.xdmf bulk_element_ids bulk_element_ids 1e-7 1e-13
-    square_5x5_tris_32_right_boundary.xdmf square_5x5_tris_32_right_boundary.xdmf bulk_node_ids bulk_node_ids 1e-7 1e-13
-    square_5x5_tris_32_left_boundary.xdmf square_5x5_tris_32_left_boundary.xdmf pressure pressure 1e-7 1e-13
-    square_5x5_tris_32_left_boundary.xdmf square_5x5_tris_32_left_boundary.xdmf bulk_element_ids bulk_element_ids 1e-7 1e-13
-    square_5x5_tris_32_left_boundary.xdmf square_5x5_tris_32_left_boundary.xdmf bulk_node_ids bulk_node_ids 1e-7 1e-13
+    square_5x5_tris_32.xdmf square_5x5_tris_32_square_5x5_tris_32.xdmf pressure pressure 1e-7 1e-13
+    square_5x5_tris_32.xdmf square_5x5_tris_32_square_5x5_tris_32.xdmf HydraulicFlow HydraulicFlow 1e-7 1e-13
+    square_5x5_tris_32.xdmf square_5x5_tris_32_square_5x5_tris_32.xdmf MaterialIDs MaterialIDs 1e-7 1e-13
+    square_5x5_tris_32.xdmf square_5x5_tris_32_square_5x5_tris_32.xdmf v v 1e-7 1e-13
+    square_5x5_tris_32_right_boundary.xdmf square_5x5_tris_32_square_5x5_tris_32_right_boundary.xdmf pressure pressure 1e-7 1e-13
+    square_5x5_tris_32_right_boundary.xdmf square_5x5_tris_32_square_5x5_tris_32_right_boundary.xdmf bulk_element_ids bulk_element_ids 1e-7 1e-13
+    square_5x5_tris_32_right_boundary.xdmf square_5x5_tris_32_square_5x5_tris_32_right_boundary.xdmf bulk_node_ids bulk_node_ids 1e-7 1e-13
+    square_5x5_tris_32_left_boundary.xdmf square_5x5_tris_32_square_5x5_tris_32_left_boundary.xdmf pressure pressure 1e-7 1e-13
+    square_5x5_tris_32_left_boundary.xdmf square_5x5_tris_32_square_5x5_tris_32_left_boundary.xdmf bulk_element_ids bulk_element_ids 1e-7 1e-13
+    square_5x5_tris_32_left_boundary.xdmf square_5x5_tris_32_square_5x5_tris_32_left_boundary.xdmf bulk_node_ids bulk_node_ids 1e-7 1e-13
 )
 
 
diff --git a/ProcessLib/Output/CreateOutput.cpp b/ProcessLib/Output/CreateOutput.cpp
index c93ad98f6baa5a096f8bb92fa6d80d1d427c2436..18df762c16118713cf469bd452ecd619ef3dd4d5 100644
--- a/ProcessLib/Output/CreateOutput.cpp
+++ b/ProcessLib/Output/CreateOutput.cpp
@@ -28,7 +28,8 @@ std::unique_ptr<Output> createOutput(
     std::vector<std::unique_ptr<MeshLib::Mesh>> const& meshes)
 {
     DBUG("Parse output configuration:");
-    OutputType const output_type = [](auto output_type) {
+    OutputType const output_type = [](auto output_type)
+    {
         try
         {
             const std::map<std::string, OutputType> outputType_to_enum = {
@@ -60,6 +61,23 @@ std::unique_ptr<Output> createOutput(
         //! \ogs_file_param{prj__time_loop__output__compress_output}
         config.getConfigParameter("compress_output", true);
 
+    auto const hdf =
+        //! \ogs_file_param{prj__time_loop__output__hdf}
+        config.getConfigSubtreeOptional("hdf");
+
+    auto number_of_files = [&hdf]() -> unsigned int
+    {
+        if (hdf)
+        {
+            //! \ogs_file_param{prj__time_loop__output__hdf__number_of_files}
+            return hdf->getConfigParameter<unsigned int>("number_of_files");
+        }
+        else
+        {
+            return 1;
+        }
+    }();
+
     auto const data_mode =
         //! \ogs_file_param{prj__time_loop__output__data_mode}
         config.getConfigParameter<std::string>("data_mode", "Appended");
@@ -156,9 +174,10 @@ std::unique_ptr<Output> createOutput(
 
     return std::make_unique<Output>(
         output_directory, output_type, prefix, suffix, compress_output,
-        data_mode, output_iteration_results, std::move(repeats_each_steps),
-        std::move(fixed_output_times), std::move(output_data_specification),
-        std::move(mesh_names_for_output), meshes);
+        number_of_files, data_mode, output_iteration_results,
+        std::move(repeats_each_steps), std::move(fixed_output_times),
+        std::move(output_data_specification), std::move(mesh_names_for_output),
+        meshes);
 }
 
 }  // namespace ProcessLib
diff --git a/ProcessLib/Output/Output.cpp b/ProcessLib/Output/Output.cpp
index c9376decfe85ea6dda33c9462b003300376653d6..6a6f80c6f6f548dbcf8211789d69e5cdb90d6cd1 100644
--- a/ProcessLib/Output/Output.cpp
+++ b/ProcessLib/Output/Output.cpp
@@ -95,7 +95,8 @@ bool Output::shallDoOutput(int timestep, double const t)
 
 Output::Output(std::string directory, OutputType file_type,
                std::string file_prefix, std::string file_suffix,
-               bool const compress_output, std::string const& data_mode,
+               bool const compress_output, unsigned int const num_of_files,
+               std::string const& data_mode,
                bool const output_nonlinear_iteration_results,
                std::vector<PairRepeatEachSteps> repeats_each_steps,
                std::vector<double>&& fixed_output_times,
@@ -107,6 +108,7 @@ Output::Output(std::string directory, OutputType file_type,
       _output_file_prefix(std::move(file_prefix)),
       _output_file_suffix(std::move(file_suffix)),
       _output_file_compression(compress_output),
+      _num_of_files(num_of_files),
       _output_file_data_mode(convertVtkDataMode(data_mode)),
       _output_nonlinear_iteration_results(output_nonlinear_iteration_results),
       _repeats_each_steps(std::move(repeats_each_steps)),
@@ -180,14 +182,16 @@ struct Output::OutputFile
                std::string const& mesh_name, int const timestep, double const t,
                int const iteration, int const data_mode_,
                bool const compression_,
-               std::set<std::string> const& outputnames)
+               std::set<std::string> const& outputnames,
+               unsigned int const num_of_files)
         : name(constructFilename(type, prefix, suffix, mesh_name, timestep, t,
                                  iteration)),
           path(BaseLib::joinPaths(directory, name)),
           type(type),
           data_mode(data_mode_),
           compression(compression_),
-          outputnames(outputnames)
+          outputnames(outputnames),
+          num_of_files(num_of_files)
     {
     }
 
@@ -199,12 +203,10 @@ struct Output::OutputFile
     /// documentation
     /// http://www.vtk.org/doc/nightly/html/classvtkXMLWriter.html
     int const data_mode;
-
     //! Enables or disables zlib-compression of the output files.
     bool const compression;
-
     std::set<std::string> outputnames;
-
+    unsigned int num_of_files;
     static std::string constructFilename(OutputType const type,
                                          std::string prefix, std::string suffix,
                                          std::string mesh_name,
@@ -248,7 +250,7 @@ void Output::outputMeshXdmf(
         _mesh_xdmf_hdf_writer = std::make_unique<MeshLib::IO::XdmfHdfWriter>(
             std::move(meshes), path, timestep, t,
             _output_data_specification.output_variables,
-            output_file.compression);
+            output_file.compression, output_file.num_of_files);
     }
     else
     {
@@ -321,7 +323,7 @@ void Output::doOutputAlways(Process const& process,
                     _output_directory, _output_file_type, _output_file_prefix,
                     _output_file_suffix, mesh.get().getName(), timestep, t,
                     iteration, _output_file_data_mode, _output_file_compression,
-                    _output_data_specification.output_variables);
+                    _output_data_specification.output_variables, 1);
 
                 pvd_file =
                     findPVDFile(process, process_id, mesh.get().getName());
@@ -331,11 +333,11 @@ void Output::doOutputAlways(Process const& process,
         else if (_output_file_type == ProcessLib::OutputType::xdmf)
         {
             std::string name = meshes[0].get().getName();
-            OutputFile const file(_output_directory, _output_file_type,
-                                  _output_file_prefix, "", name, timestep, t,
-                                  iteration, _output_file_data_mode,
-                                  _output_file_compression,
-                                  _output_data_specification.output_variables);
+            OutputFile const file(
+                _output_directory, _output_file_type, _output_file_prefix, "",
+                name, timestep, t, iteration, _output_file_data_mode,
+                _output_file_compression,
+                _output_data_specification.output_variables, _num_of_files);
 
             outputMeshXdmf(std::move(file), std::move(meshes), timestep, t);
         }
diff --git a/ProcessLib/Output/Output.h b/ProcessLib/Output/Output.h
index 6e530c7f220dbe3e33034f5f747ea905bbb24eb6..075f9f98241eb3470a7252ef8147c0f805e80ae3 100644
--- a/ProcessLib/Output/Output.h
+++ b/ProcessLib/Output/Output.h
@@ -40,7 +40,7 @@ public:
 public:
     Output(std::string directory, OutputType const type, std::string prefix,
            std::string suffix, bool const compress_output,
-           std::string const& data_mode,
+           unsigned int num_of_files, std::string const& data_mode,
            bool const output_nonlinear_iteration_results,
            std::vector<PairRepeatEachSteps> repeats_each_steps,
            std::vector<double>&& fixed_output_times,
@@ -108,6 +108,8 @@ private:
 
     //! Enables or disables zlib-compression of the output files.
     bool const _output_file_compression;
+    //! Specifies the number of hdf5 output files.
+    unsigned int const _num_of_files;
 
     //! Chooses vtk's data mode for output following the enumeration given in
     /// the vtkXMLWriter: {Ascii, Binary, Appended}.  See vtkXMLWriter
diff --git a/ProcessLib/SteadyStateDiffusion/Tests.cmake b/ProcessLib/SteadyStateDiffusion/Tests.cmake
index 688763e93af237bb9a713c0c1f4e1da8d851e9ee..8ca8c99f3c7a5f4bd818bbf76190fbf32d19b99f 100644
--- a/ProcessLib/SteadyStateDiffusion/Tests.cmake
+++ b/ProcessLib/SteadyStateDiffusion/Tests.cmake
@@ -425,8 +425,36 @@ AddTest(
     TESTER xdmfdiff
     REQUIREMENTS OGS_USE_MPI
     DIFF_DATA
-    cube_1e3_np3.xdmf cube_1x1x1_hex_1e3.xdmf pressure pressure 1e-3 1e-3
-    cube_1e3_np3.xdmf cube_1x1x1_hex_1e3.xdmf v v 1e-3 1e-3
+    cube_1e3_np3.xdmf cube_1e3_np3_cube_1x1x1_hex_1e3.xdmf pressure pressure 1e-3 1e-3
+    cube_1e3_np3.xdmf cube_1e3_np3_cube_1x1x1_hex_1e3.xdmf v v 1e-3 1e-3
+)
+
+AddTest(
+    NAME ParallelFEM_GroundWaterFlow3D_NeumannBC_XDMF_np3_2files
+    PATH EllipticPETSc
+    EXECUTABLE ogs
+    EXECUTABLE_ARGS cube_1e3_XDMF_np3_2files.prj
+    WRAPPER mpirun
+    WRAPPER_ARGS -np 3
+    TESTER xdmfdiff
+    REQUIREMENTS OGS_USE_MPI
+    DIFF_DATA
+    cube_1e3_np3_2files_0.xdmf cube_1e3_np3_cube_1x1x1_hex_1e3.xdmf pressure pressure 1e-3 1e-3
+    cube_1e3_np3_2files_0.xdmf cube_1e3_np3_cube_1x1x1_hex_1e3.xdmf v v 1e-3 1e-3
+)
+
+AddTest(
+    NAME ParallelFEM_GroundWaterFlow3D_NeumannBC_XDMF_np3_3files
+    PATH EllipticPETSc
+    EXECUTABLE ogs
+    EXECUTABLE_ARGS cube_1e3_XDMF_np3_3files.prj
+    WRAPPER mpirun
+    WRAPPER_ARGS -np 3
+    TESTER xdmfdiff
+    REQUIREMENTS OGS_USE_MPI
+    DIFF_DATA
+    cube_1e3_np3_3files_0.xdmf cube_1e3_np3_cube_1x1x1_hex_1e3.xdmf pressure pressure 1e-3 1e-3
+    cube_1e3_np3_3files_0.xdmf cube_1e3_np3_cube_1x1x1_hex_1e3.xdmf v v 1e-3 1e-3
 )
 
 AddTest(
@@ -439,8 +467,8 @@ AddTest(
     TESTER xdmfdiff
     REQUIREMENTS OGS_USE_MPI
     DIFF_DATA
-    cube_1e3_np2.xdmf cube_1x1x1_hex_1e3.xdmf pressure pressure 1e-3 1e-3
-    cube_1e3_np2.xdmf cube_1x1x1_hex_1e3.xdmf v v 1e-3 1e-3
+    cube_1e3_np2.xdmf cube_1e3_np2_cube_1x1x1_hex_1e3.xdmf pressure pressure 1e-3 1e-3
+    cube_1e3_np2.xdmf cube_1e3_np2_cube_1x1x1_hex_1e3.xdmf v v 1e-3 1e-3
 )
 
 AddTest(
diff --git a/Tests/Data/EllipticPETSc/cube_1e3_XDMF_np3.prj b/Tests/Data/EllipticPETSc/cube_1e3_XDMF_np3.prj
index 74c0442b221169a2add593fbc652967a4940882d..3fb3aef3a66e5fd5b1cdf8cd8e95de8d567f7714 100644
--- a/Tests/Data/EllipticPETSc/cube_1e3_XDMF_np3.prj
+++ b/Tests/Data/EllipticPETSc/cube_1e3_XDMF_np3.prj
@@ -65,6 +65,9 @@
                 <variable> v      </variable>
             </variables>
             <suffix>_ts_{:timestep}_t_{:time}</suffix>
+            <hdf>
+                <num_of_files>1</num_of_files>
+            </hdf>
         </output>
     </time_loop>
     <parameters>