diff --git a/CMakeLists.txt b/CMakeLists.txt
index cc718b9c0476e0726dacde7e6a7e0939b861b368..89463224a89678ee16c84af00e469d88ec185054 100644
--- a/CMakeLists.txt
+++ b/CMakeLists.txt
@@ -137,4 +137,5 @@ message("    MPI:     ${ADIOS_USE_MPI}")
 message("    BZip2:   ${ADIOS_USE_BZip2}")
 message("    ADIOS1:  ${ADIOS_USE_ADIOS1}")
 message("    DataMan: ${ADIOS_USE_DataMan}")
+message("    HDF5: ${ADIOS_USE_PHDF5}")
 message("")
diff --git a/examples/hello/CMakeLists.txt b/examples/hello/CMakeLists.txt
index 0cf0a5df520d8388fa7c5e19f9ca8c9a94aec614..73fa0f560c36f0c1c0f1744cb44a90f7f7df6ba2 100644
--- a/examples/hello/CMakeLists.txt
+++ b/examples/hello/CMakeLists.txt
@@ -14,3 +14,7 @@ if(ADIOS_USE_DataMan)
   add_subdirectory(datamanReader)
   add_subdirectory(datamanWriter)
 endif()
+
+if(ADIOS_USE_PHDF5)
+  add_subdirectory(hdf5Writer)
+endif()
diff --git a/examples/hello/hdf5Writer/CMakeLists.txt b/examples/hello/hdf5Writer/CMakeLists.txt
new file mode 100644
index 0000000000000000000000000000000000000000..2f7d7fc0cbc69fc167f32a60b4aa030b44174bd1
--- /dev/null
+++ b/examples/hello/hdf5Writer/CMakeLists.txt
@@ -0,0 +1,30 @@
+#------------------------------------------------------------------------------#
+# Distributed under the OSI-approved Apache License, Version 2.0.  See
+# accompanying file Copyright.txt for details.
+#------------------------------------------------------------------------------#
+
+#add_executable(hello_hdf5Writer_nompi helloHDF5Writer_nompi.cpp)
+message("    HDF5:   ${HDF5_INCLUDE_DIRS}")
+include_directories(${HDF5_INCLUDE_DIRS})
+#target_link_libraries(hello_hdf5Writer_nompi adios2_nompi)
+
+if(ADIOS_BUILD_TESTING)
+  #add_test(
+  #  NAME Example::hello::hdf5Writer_nompi
+  #  COMMAND hello_hdf5Writer_nompi
+  #)
+endif()
+
+if(ADIOS_USE_MPI)
+  find_package(MPI COMPONENTS C REQUIRED)
+
+  add_executable(hello_hdf5Writer helloHDF5Writer.cpp)
+  target_link_libraries(hello_hdf5Writer adios2)
+  target_include_directories(hello_hdf5Writer PRIVATE ${MPI_C_INCLUDE_PATH})
+  target_link_libraries(hello_hdf5Writer adios2 ${MPI_C_LIBRARIES})
+
+  if(ADIOS_BUILD_TESTING)
+    add_test(NAME Example::hello::hdf5Writer COMMAND hello_hdf5Writer)
+  endif()
+endif()
+
diff --git a/examples/hello/hdf5Writer/helloHDF5Writer.cpp b/examples/hello/hdf5Writer/helloHDF5Writer.cpp
new file mode 100644
index 0000000000000000000000000000000000000000..3b32828404b86a4e3d8b7bda77a33cdce16a141b
--- /dev/null
+++ b/examples/hello/hdf5Writer/helloHDF5Writer.cpp
@@ -0,0 +1,144 @@
+/*
+ * HDF5Writer.cpp
+ *
+ *  Created on: March 20, 2017
+ *      Author: Junmin
+ */
+
+#include <iostream>
+#include <vector>
+
+#include <mpi.h>
+
+#include "ADIOS_CPP.h"
+
+int main(int argc, char *argv[])
+{
+    MPI_Init(&argc, &argv);
+    int rank, size;
+    MPI_Comm_rank(MPI_COMM_WORLD, &rank);
+    MPI_Comm_size(MPI_COMM_WORLD, &size);
+
+    const bool adiosDebug = true;
+    adios::ADIOS adios(MPI_COMM_WORLD, adios::Verbose::INFO, adiosDebug);
+
+    // Application variable
+    const std::size_t intDim1 = 4;
+    const std::size_t intDim2 = 3;
+    std::vector<int> myInts = {10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, 21};
+
+    std::vector<double> myDoubles = {0, 1, 2, 3, 4, 5, 6, 7, 8, 9};
+    const std::size_t Nx = myDoubles.size();
+
+    std::vector<std::complex<float>> myCFloats;
+    const std::size_t ComplexDataSize = 3;
+    myCFloats.reserve(ComplexDataSize);
+    myCFloats.emplace_back(1, 3);
+    myCFloats.emplace_back(2, 2);
+    myCFloats.emplace_back(3, 1);
+
+    std::vector<std::complex<double>> myCDoubles;
+    myCDoubles.reserve(ComplexDataSize);
+    myCDoubles.emplace_back(1.1, -3.3);
+    myCDoubles.emplace_back(2.1, -2.2);
+    myCDoubles.emplace_back(3.1, -1.1);
+
+    std::vector<std::complex<long double>> myCLongDoubles;
+    myCLongDoubles.reserve(ComplexDataSize);
+    myCLongDoubles.emplace_back(1.11, -3.33);
+    myCLongDoubles.emplace_back(2.11, -2.22);
+    myCLongDoubles.emplace_back(3.11, -1.11);
+
+    std::size_t doubleVCount = Nx / size;
+    std::size_t complexCount = ComplexDataSize / size;
+    std::size_t intCountDim1 = intDim1 / size;
+
+    std::size_t doubleVOffset = rank * doubleVCount;
+    std::size_t complexOffset = rank * complexCount;
+    std::size_t intOffsetDim1 = rank * intCountDim1;
+    std::size_t intOffsetDim2 = 0;
+
+    if ((size > 1) && (rank == size - 1))
+    {
+        doubleVCount = Nx - rank * (Nx / size);
+        complexCount = ComplexDataSize - rank * (ComplexDataSize / size);
+        intCountDim1 = intDim1 - rank * (intDim1 / size);
+    }
+
+    try
+    {
+        // Define variable and local size
+        auto &ioMyInts =
+            adios.DefineVariable<int>("myInts", {intCountDim1, intDim2}, {4, 3},
+                                      {intOffsetDim1, intOffsetDim2});
+        auto &ioMyDoubles = adios.DefineVariable<double>(
+            "myDoubles", {doubleVCount}, {Nx}, {doubleVOffset});
+        auto &ioMyCFloats = adios.DefineVariable<std::complex<float>>(
+            "myCFloats", {complexCount}, {3}, {complexOffset});
+        auto &ioMyCDoubles = adios.DefineVariable<std::complex<double>>(
+            "myCDoubles", {complexCount}, {3}, {complexOffset});
+        auto &ioMyCLongDoubles =
+            adios.DefineVariable<std::complex<long double>>(
+                "myCLongDoubles", {complexCount}, {3}, {complexOffset});
+
+        // Define method for engine creation, it is basically straight-forward
+        // parameters
+        adios::Method &HDF5Settings = adios.DeclareMethod("w");
+        HDF5Settings.SetEngine("HDF5Writer");
+        HDF5Settings.SetParameters("chunck=yes", "collectiveIO=yes");
+        // HDF5Settings.AddTransport( "Mdtm", "localIP=128.0.0.0.1",
+        // "remoteIP=128.0.0.0.2", "tolerances=1,2,3" );
+
+        // Create engine smart pointer to HDF5 Engine due to polymorphism,
+        // Open returns a smart pointer to Engine containing the Derived class
+        // HDF5
+        auto HDF5Writer = adios.Open("test.h5", "w", HDF5Settings);
+
+        if (HDF5Writer == nullptr)
+            throw std::ios_base::failure(
+                "ERROR: failed to create HDF5 I/O engine at Open\n");
+
+        HDF5Writer->Write(ioMyDoubles, myDoubles.data() +
+                                           doubleVOffset); // Base class Engine
+                                                           // own the Write<T>
+                                                           // that will call
+                                                           // overloaded Write
+                                                           // from Derived
+        HDF5Writer->Write(ioMyInts,
+                          myInts.data() + (intOffsetDim1 * intDim2 * rank));
+
+        HDF5Writer->Write(ioMyCFloats, myCFloats.data() + complexOffset);
+        HDF5Writer->Write(ioMyCDoubles, myCDoubles.data() + complexOffset);
+        HDF5Writer->Write(ioMyCLongDoubles,
+                          myCLongDoubles.data() + complexOffset);
+        HDF5Writer->Close();
+    }
+    catch (std::invalid_argument &e)
+    {
+        if (rank == 0)
+        {
+            std::cout << "Invalid argument exception, STOPPING PROGRAM\n";
+            std::cout << e.what() << "\n";
+        }
+    }
+    catch (std::ios_base::failure &e)
+    {
+        if (rank == 0)
+        {
+            std::cout << "System exception, STOPPING PROGRAM\n";
+            std::cout << e.what() << "\n";
+        }
+    }
+    catch (std::exception &e)
+    {
+        if (rank == 0)
+        {
+            std::cout << "Exception, STOPPING PROGRAM\n";
+            std::cout << e.what() << "\n";
+        }
+    }
+
+    MPI_Finalize();
+
+    return 0;
+}
diff --git a/include/engine/hdf5/HDF5ReaderP.h b/include/engine/hdf5/HDF5ReaderP.h
new file mode 100644
index 0000000000000000000000000000000000000000..ec5d910ebf453c4f4f9eb6d5faf55248759feef6
--- /dev/null
+++ b/include/engine/hdf5/HDF5ReaderP.h
@@ -0,0 +1,14 @@
+/*
+ * Distributed under the OSI-approved Apache License, Version 2.0.  See
+ * accompanying file Copyright.txt for details.
+ *
+ * HDF5ReaderP.h
+ *
+ *  Created on: March 20, 2017
+ *      Author: Junmin
+ */
+
+#ifndef HDF5_READER_P_H
+#define HDF5_READER_P_H
+
+#endif
diff --git a/include/engine/hdf5/HDF5WriterP.h b/include/engine/hdf5/HDF5WriterP.h
new file mode 100644
index 0000000000000000000000000000000000000000..cbdbf1c7dcdd5358a11386a29c9bff1f58f606fd
--- /dev/null
+++ b/include/engine/hdf5/HDF5WriterP.h
@@ -0,0 +1,113 @@
+
+/*
+ * Distributed under the OSI-approved Apache License, Version 2.0.  See
+ * accompanying file Copyright.txt for details.
+ *
+ * HDF5WriterP.h
+ *
+ *  Created on: March 20, 2017
+ *      Author: Junmin
+ */
+
+#ifndef HDF5_WRITER_P_H_
+#define HDF5_WRITER_P_H_
+
+#include "core/Engine.h"
+
+// supported capsules
+#include "capsule/heap/STLVector.h"
+
+#include "ADIOS_MPI.h"
+
+#include <hdf5.h>
+
+namespace adios
+{
+
+class HDF5Writer : public Engine
+{
+
+public:
+    /**
+     * Constructor for single BP capsule engine, writes in BP format into a
+     * single
+     * heap capsule
+     * @param name unique name given to the engine
+     * @param accessMode
+     * @param mpiComm
+     * @param method
+     * @param debugMode
+     */
+    HDF5Writer(ADIOS &adios, const std::string name,
+               const std::string accessMode, MPI_Comm mpiComm,
+               const Method &method);
+
+    virtual ~HDF5Writer();
+
+    void Write(Variable<char> &variable, const char *values);
+    void Write(Variable<unsigned char> &variable, const unsigned char *values);
+    void Write(Variable<short> &variable, const short *values);
+    void Write(Variable<unsigned short> &variable,
+               const unsigned short *values);
+    void Write(Variable<int> &variable, const int *values);
+    void Write(Variable<unsigned int> &variable, const unsigned int *values);
+    void Write(Variable<long int> &variable, const long int *values);
+    void Write(Variable<unsigned long int> &variable,
+               const unsigned long int *values);
+    void Write(Variable<long long int> &variable, const long long int *values);
+    void Write(Variable<unsigned long long int> &variable,
+               const unsigned long long int *values);
+    void Write(Variable<float> &variable, const float *values);
+    void Write(Variable<double> &variable, const double *values);
+    void Write(Variable<long double> &variable, const long double *values);
+    void Write(Variable<std::complex<float>> &variable,
+               const std::complex<float> *values);
+    void Write(Variable<std::complex<double>> &variable,
+               const std::complex<double> *values);
+    void Write(Variable<std::complex<long double>> &variable,
+               const std::complex<long double> *values);
+
+    void Write(const std::string variableName, const char *values);
+    void Write(const std::string variableName, const unsigned char *values);
+    void Write(const std::string variableName, const short *values);
+    void Write(const std::string variableName, const unsigned short *values);
+    void Write(const std::string variableName, const int *values);
+    void Write(const std::string variableName, const unsigned int *values);
+    void Write(const std::string variableName, const long int *values);
+    void Write(const std::string variableName, const unsigned long int *values);
+    void Write(const std::string variableName, const long long int *values);
+    void Write(const std::string variableName,
+               const unsigned long long int *values);
+    void Write(const std::string variableName, const float *values);
+    void Write(const std::string variableName, const double *values);
+    void Write(const std::string variableName, const long double *values);
+    void Write(const std::string variableName,
+               const std::complex<float> *values);
+    void Write(const std::string variableName,
+               const std::complex<double> *values);
+    void Write(const std::string variableName,
+               const std::complex<long double> *values);
+
+    void Close(const int transportIndex = -1);
+
+private:
+    capsule::STLVector
+        m_Buffer; ///< heap capsule, contains data and metadata buffers
+
+    void Init();
+    void clean();
+
+    hid_t _plist_id, _file_id, _dset_id;
+    hid_t _memspace, _filespace;
+
+    hid_t DefH5T_COMPLEX_DOUBLE;
+    hid_t DefH5T_COMPLEX_FLOAT;
+    hid_t DefH5T_COMPLEX_LongDOUBLE;
+
+    template <class T>
+    void UseHDFWrite(Variable<T> &variable, const T *values, hid_t h5type);
+};
+
+} // end namespace adios
+
+#endif /* HDF5_WRITER_P_H_ */
diff --git a/source/ADIOS.cpp b/source/ADIOS.cpp
index 2e571a48f1a9b95ee66b2f031973c5c6db77ba32..30aac797cd4d6533626106baef0dcab3eeed529b 100644
--- a/source/ADIOS.cpp
+++ b/source/ADIOS.cpp
@@ -36,6 +36,12 @@
 #include "engine/adios1/ADIOS1Writer.h"
 #endif
 
+#ifdef ADIOS_HAVE_PHDF5 // external dependencies
+#ifdef ADIOS_HAVE_MPI
+#include "engine/hdf5/HDF5ReaderP.h"
+#include "engine/hdf5/HDF5WriterP.h"
+#endif
+#endif
 namespace adios
 {
 
@@ -186,6 +192,17 @@ std::shared_ptr<Engine> ADIOS::Open(const std::string &name,
         // method,
         // iomode, timeout_sec, m_DebugMode, method.m_nThreads );
     }
+    else if (type == "HDF5Writer") // -junmin
+    {
+#if defined(ADIOS_HAVE_PHDF5) && defined(ADIOS_HAVE_MPI)
+        return std::make_shared<HDF5Writer>(*this, name, accessMode, mpiComm,
+                                            method);
+#else
+        throw std::invalid_argument("ERROR: this version didn't compile with "
+                                    "HDF5 library, can't use HDF5\n");
+#endif
+    }
+
     else
     {
         if (m_DebugMode == true)
diff --git a/source/CMakeLists.txt b/source/CMakeLists.txt
index b01125defc82fcde37754d67f12756fd530875ad..62a642e966648cd5d2d73c4473146429920fa58a 100644
--- a/source/CMakeLists.txt
+++ b/source/CMakeLists.txt
@@ -25,6 +25,12 @@ foreach(adios2_target IN LISTS adios2_targets)
   
     engine/bp/BPFileReader.cpp
     engine/bp/BPFileWriter.cpp
+
+    utilities/format/bp1/BP1Base.cpp
+    utilities/format/bp1/BP1Aggregator.cpp
+    utilities/format/bp1/BP1Writer.cpp
+    
+    utilities/profiling/iochrono/Timer.cpp
   
     functions/adiosFunctions.cpp
   
@@ -54,6 +60,7 @@ foreach(adios2_target IN LISTS adios2_targets)
     target_compile_definitions(${adios2_target} PRIVATE ADIOS_HAVE_DATAMAN)
     target_link_libraries(${adios2_target} PRIVATE ${CMAKE_DL_LIBS})
   endif()
+
   
   if(ADIOS_USE_BZip2)
     find_package(BZip2 REQUIRED)
@@ -83,4 +90,20 @@ if(ADIOS_USE_MPI)
     target_compile_definitions(adios2 PRIVATE ADIOS_HAVE_ADIOS1)
     target_link_libraries(adios2 PRIVATE adios1::adios)
   endif()
+
+  if(ADIOS_USE_PHDF5)
+    find_package(HDF5 REQUIRED)
+    message("    HDF5 ROOT:   ${HDF5_ROOT}")
+    message("    HDF5 include:   ${HDF5_INCLUDE_DIRS}   is paralle? ${HDF5_IS_PARALLEL}")
+
+    target_include_directories(adios2 PRIVATE ${HDF5_INCLUDE_DIRS})
+
+    target_sources(adios2 PRIVATE
+      engine/hdf5/HDF5ReaderP.cpp
+      engine/hdf5/HDF5WriterP.cpp
+    )
+    target_compile_definitions(adios2 PRIVATE ADIOS_HAVE_PHDF5)
+    target_link_libraries(adios2 PRIVATE ${HDF5_C_LIBRARIES})
+  endif()
+
 endif()
diff --git a/source/engine/hdf5/HDF5ReaderP.cpp b/source/engine/hdf5/HDF5ReaderP.cpp
new file mode 100644
index 0000000000000000000000000000000000000000..00ec0a8fbc3098e3c12c8647eb3113d4bd5d69f0
--- /dev/null
+++ b/source/engine/hdf5/HDF5ReaderP.cpp
@@ -0,0 +1,13 @@
+/*
+ * Distributed under the OSI-approved Apache License, Version 2.0.  See
+ * accompanying file Copyright.txt for details.
+ *
+ * HDF5WriterP.cpp
+ *
+ *  Created on: March 20, 2017
+ *      Author: Junmin
+ */
+
+#include <iostream> //needs to go away, this is just for demo purposes
+
+#include "engine/hdf5/HDF5ReaderP.h"
diff --git a/source/engine/hdf5/HDF5WriterP.cpp b/source/engine/hdf5/HDF5WriterP.cpp
new file mode 100644
index 0000000000000000000000000000000000000000..05656e2babf1170f368702b2cd0d309d47800a8c
--- /dev/null
+++ b/source/engine/hdf5/HDF5WriterP.cpp
@@ -0,0 +1,365 @@
+/*
+ * Distributed under the OSI-approved Apache License, Version 2.0.  See
+ * accompanying file Copyright.txt for details.
+ *
+ * HDF5WriterP.cpp
+ *
+ *  Created on: March 20, 2017
+ *      Author: Junmin
+ */
+
+#include <iostream> //needs to go away, this is just for demo purposes
+
+#include "engine/hdf5/HDF5WriterP.h"
+
+#include "core/Support.h"
+#include "functions/adiosFunctions.h" //CSVToVector
+
+namespace adios
+{
+
+HDF5Writer::HDF5Writer(ADIOS &adios, const std::string name,
+                       const std::string accessMode, MPI_Comm mpiComm,
+                       const Method &method)
+: Engine(adios, "HDF5Writer", name, accessMode, mpiComm,
+         method, /*debugMode, cores,*/
+         " HDF5Writer constructor (or call to ADIOS Open).\n"),
+  m_Buffer(accessMode, m_RankMPI, m_DebugMode)
+{
+    //
+    //  16, 4 vs: 8
+    // std::cout<<sizeof(std::complex<double>)<<",
+    // "<<sizeof(H5T_NATIVE_DOUBLE)<<" vs:
+    // "<<H5Tget_size(H5T_NATIVE_DOUBLE)<<std::endl;
+    //  8, 4 vs: 4
+    // std::cout<<sizeof(std::complex<float>)<<", "<<sizeof(H5T_NATIVE_FLOAT)<<"
+    // vs: "<<H5Tget_size(H5T_NATIVE_FLOAT)<<std::endl;
+    //  32, 4 vs: 16
+    // std::cout<<sizeof(std::complex<long double>)<<",
+    // "<<sizeof(H5T_NATIVE_LDOUBLE)<<" vs:
+    // "<<H5Tget_size(H5T_NATIVE_LDOUBLE)<<std::endl;
+
+    DefH5T_COMPLEX_FLOAT = H5Tcreate(H5T_COMPOUND, sizeof(std::complex<float>));
+    H5Tinsert(DefH5T_COMPLEX_FLOAT, "freal", 0, H5T_NATIVE_FLOAT);
+    H5Tinsert(DefH5T_COMPLEX_FLOAT, "fimg", H5Tget_size(H5T_NATIVE_FLOAT),
+              H5T_NATIVE_FLOAT);
+
+    DefH5T_COMPLEX_DOUBLE =
+        H5Tcreate(H5T_COMPOUND, sizeof(std::complex<double>));
+    H5Tinsert(DefH5T_COMPLEX_DOUBLE, "dreal", 0, H5T_NATIVE_DOUBLE);
+    H5Tinsert(DefH5T_COMPLEX_DOUBLE, "dimg", H5Tget_size(H5T_NATIVE_DOUBLE),
+              H5T_NATIVE_DOUBLE);
+
+    DefH5T_COMPLEX_LongDOUBLE =
+        H5Tcreate(H5T_COMPOUND, sizeof(std::complex<long double>));
+    H5Tinsert(DefH5T_COMPLEX_LongDOUBLE, "ldouble real", 0, H5T_NATIVE_LDOUBLE);
+    H5Tinsert(DefH5T_COMPLEX_LongDOUBLE, "ldouble img",
+              H5Tget_size(H5T_NATIVE_LDOUBLE), H5T_NATIVE_LDOUBLE);
+
+    Init();
+}
+
+HDF5Writer::~HDF5Writer() {}
+
+void HDF5Writer::Init()
+{
+    if (m_AccessMode != "w" && m_AccessMode != "write" && m_AccessMode != "a" &&
+        m_AccessMode != "append")
+    {
+        throw std::invalid_argument(
+            "ERROR: HDF5Writer doesn't support access mode " + m_AccessMode +
+            ", in call to ADIOS Open or HDF5Writer constructor\n");
+    }
+    // std::cout << "method: # of inputs:" << m_Method.m_Parameters.size() <<
+    // std::endl;
+
+    // std::cout << "::Init hdf5 parallel writer. File name:" << m_Name <<
+    // std::endl;
+
+    _plist_id = H5Pcreate(H5P_FILE_ACCESS);
+
+    H5Pset_fapl_mpio(_plist_id, m_MPIComm, MPI_INFO_NULL);
+
+    /*
+     * Create a new file collectively and release property list identifier.
+     */
+    _file_id = H5Fcreate(m_Name.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, _plist_id);
+    H5Pclose(_plist_id);
+}
+
+void HDF5Writer::Write(Variable<char> &variable, const char *values)
+{
+    UseHDFWrite(variable, values, H5T_NATIVE_CHAR);
+}
+
+void HDF5Writer::Write(Variable<unsigned char> &variable,
+                       const unsigned char *values)
+{
+    UseHDFWrite(variable, values, H5T_NATIVE_UCHAR);
+}
+
+void HDF5Writer::Write(Variable<short> &variable, const short *values)
+{
+    UseHDFWrite(variable, values, H5T_NATIVE_SHORT);
+}
+
+void HDF5Writer::Write(Variable<unsigned short> &variable,
+                       const unsigned short *values)
+{
+    UseHDFWrite(variable, values, H5T_NATIVE_USHORT);
+}
+
+void HDF5Writer::Write(Variable<int> &variable, const int *values)
+{
+    UseHDFWrite(variable, values, H5T_NATIVE_INT);
+}
+
+void HDF5Writer::Write(Variable<unsigned int> &variable,
+                       const unsigned int *values)
+{
+    UseHDFWrite(variable, values, H5T_NATIVE_UINT);
+}
+
+void HDF5Writer::Write(Variable<long int> &variable, const long int *values)
+{
+    UseHDFWrite(variable, values, H5T_NATIVE_LONG);
+}
+
+void HDF5Writer::Write(Variable<unsigned long int> &variable,
+                       const unsigned long int *values)
+{
+    UseHDFWrite(variable, values, H5T_NATIVE_ULONG);
+}
+
+void HDF5Writer::Write(Variable<long long int> &variable,
+                       const long long int *values)
+{
+    UseHDFWrite(variable, values, H5T_NATIVE_LLONG);
+}
+
+void HDF5Writer::Write(Variable<unsigned long long int> &variable,
+                       const unsigned long long int *values)
+{
+    UseHDFWrite(variable, values, H5T_NATIVE_ULLONG);
+}
+
+void HDF5Writer::Write(Variable<float> &variable, const float *values)
+{
+    UseHDFWrite(variable, values, H5T_NATIVE_FLOAT);
+}
+
+void HDF5Writer::Write(Variable<double> &variable, const double *values)
+{
+    UseHDFWrite(variable, values, H5T_NATIVE_DOUBLE);
+}
+
+void HDF5Writer::Write(Variable<long double> &variable,
+                       const long double *values)
+{
+    UseHDFWrite(variable, values, H5T_NATIVE_LDOUBLE);
+}
+
+void HDF5Writer::Write(Variable<std::complex<float>> &variable,
+                       const std::complex<float> *values)
+{
+    UseHDFWrite(variable, values, DefH5T_COMPLEX_FLOAT);
+}
+
+void HDF5Writer::Write(Variable<std::complex<double>> &variable,
+                       const std::complex<double> *values)
+{
+    UseHDFWrite(variable, values, DefH5T_COMPLEX_DOUBLE);
+}
+
+void HDF5Writer::Write(Variable<std::complex<long double>> &variable,
+                       const std::complex<long double> *values)
+{
+    UseHDFWrite(variable, values, DefH5T_COMPLEX_LongDOUBLE);
+}
+
+// String version
+void HDF5Writer::Write(const std::string variableName, const char *values)
+{
+    UseHDFWrite(m_ADIOS.GetVariable<char>(variableName), values,
+                H5T_NATIVE_CHAR);
+}
+
+void HDF5Writer::Write(const std::string variableName,
+                       const unsigned char *values)
+{
+    UseHDFWrite(m_ADIOS.GetVariable<unsigned char>(variableName), values,
+                H5T_NATIVE_UCHAR);
+}
+
+void HDF5Writer::Write(const std::string variableName, const short *values)
+{
+    UseHDFWrite(m_ADIOS.GetVariable<short>(variableName), values,
+                H5T_NATIVE_SHORT);
+}
+
+void HDF5Writer::Write(const std::string variableName,
+                       const unsigned short *values)
+{
+    UseHDFWrite(m_ADIOS.GetVariable<unsigned short>(variableName), values,
+                H5T_NATIVE_USHORT);
+}
+
+void HDF5Writer::Write(const std::string variableName, const int *values)
+{
+    UseHDFWrite(m_ADIOS.GetVariable<int>(variableName), values, H5T_NATIVE_INT);
+}
+
+void HDF5Writer::Write(const std::string variableName,
+                       const unsigned int *values)
+{
+    UseHDFWrite(m_ADIOS.GetVariable<unsigned int>(variableName), values,
+                H5T_NATIVE_UINT);
+}
+
+void HDF5Writer::Write(const std::string variableName, const long int *values)
+{
+    UseHDFWrite(m_ADIOS.GetVariable<long int>(variableName), values,
+                H5T_NATIVE_LONG);
+}
+
+void HDF5Writer::Write(const std::string variableName,
+                       const unsigned long int *values)
+{
+    UseHDFWrite(m_ADIOS.GetVariable<unsigned long int>(variableName), values,
+                H5T_NATIVE_ULONG);
+}
+
+void HDF5Writer::Write(const std::string variableName,
+                       const long long int *values)
+{
+    UseHDFWrite(m_ADIOS.GetVariable<long long int>(variableName), values,
+                H5T_NATIVE_LLONG);
+}
+
+void HDF5Writer::Write(const std::string variableName,
+                       const unsigned long long int *values)
+{
+    UseHDFWrite(m_ADIOS.GetVariable<unsigned long long int>(variableName),
+                values, H5T_NATIVE_ULLONG);
+}
+
+void HDF5Writer::Write(const std::string variableName, const float *values)
+{
+    UseHDFWrite(m_ADIOS.GetVariable<float>(variableName), values,
+                H5T_NATIVE_FLOAT);
+}
+
+void HDF5Writer::Write(const std::string variableName, const double *values)
+{
+    UseHDFWrite(m_ADIOS.GetVariable<double>(variableName), values,
+                H5T_NATIVE_DOUBLE);
+}
+
+void HDF5Writer::Write(const std::string variableName,
+                       const long double *values)
+{
+    UseHDFWrite(m_ADIOS.GetVariable<long double>(variableName), values,
+                H5T_NATIVE_LDOUBLE);
+}
+
+void HDF5Writer::Write(const std::string variableName,
+                       const std::complex<float> *values)
+{
+    UseHDFWrite(m_ADIOS.GetVariable<std::complex<float>>(variableName), values,
+                DefH5T_COMPLEX_FLOAT);
+}
+
+void HDF5Writer::Write(const std::string variableName,
+                       const std::complex<double> *values)
+{
+    UseHDFWrite(m_ADIOS.GetVariable<std::complex<double>>(variableName), values,
+                DefH5T_COMPLEX_DOUBLE);
+}
+
+void HDF5Writer::Write(const std::string variableName,
+                       const std::complex<long double> *values)
+{
+    UseHDFWrite(m_ADIOS.GetVariable<std::complex<long double>>(variableName),
+                values, DefH5T_COMPLEX_LongDOUBLE);
+}
+
+void HDF5Writer::Close(const int transportIndex)
+{
+    // std::cout << " ===> CLOSING HDF5 <===== " << std::endl;
+    // H5Dclose(_dset_id);
+    // H5Sclose(_filespace);
+    // H5Sclose(_memspace);
+    // H5Pclose(_plist_id);
+    H5Fclose(_file_id);
+}
+
+template <class T>
+void HDF5Writer::UseHDFWrite(Variable<T> &variable, const T *values,
+                             hid_t h5type)
+{
+    // here comes your magic at Writing now variable.m_UserValues has the data
+    // passed by the user
+    // set variable
+    variable.m_AppValues = values;
+    m_WrittenVariables.insert(variable.m_Name);
+
+    int dimSize = variable.m_GlobalDimensions.size();
+    /*
+    std::cout << "writting : " << variable.m_Name
+              << " dim size:" << variable.m_GlobalDimensions.size() <<
+    std::endl;
+
+    for (int i = 0; i < dimSize; i++) {
+      std::cout << " dim: " << i << ", size:" << variable.m_GlobalDimensions[i]
+                << " offset=" << variable.m_Offsets[i]
+                << " count=" << variable.m_LocalDimensions[i] << std::endl;
+    }
+    */
+    std::vector<hsize_t> dimsf, count, offset;
+
+    for (int i = 0; i < dimSize; i++)
+    {
+        dimsf.push_back(variable.m_GlobalDimensions[i]);
+        count.push_back(variable.m_LocalDimensions[i]);
+        offset.push_back(variable.m_Offsets[i]);
+    }
+
+    _filespace = H5Screate_simple(dimSize, dimsf.data(), NULL);
+
+    _dset_id = H5Dcreate(_file_id, variable.m_Name.c_str(), h5type, _filespace,
+                         H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
+    // H5Sclose(_filespace);
+
+    _memspace = H5Screate_simple(dimSize, count.data(), NULL);
+
+    // Select hyperslab
+    _filespace = H5Dget_space(_dset_id);
+    H5Sselect_hyperslab(_filespace, H5S_SELECT_SET, offset.data(), NULL,
+                        count.data(), NULL);
+
+    //  Create property list for collective dataset write.
+
+    _plist_id = H5Pcreate(H5P_DATASET_XFER);
+    H5Pset_dxpl_mpio(_plist_id, H5FD_MPIO_COLLECTIVE);
+
+    herr_t status;
+
+    status =
+        H5Dwrite(_dset_id, h5type, _memspace, _filespace, _plist_id, values);
+
+    if (status < 0)
+    {
+        // error
+        std::cerr << " Write failed. " << std::endl;
+    }
+
+    // std::cout << " ==> User is responsible for freeing the data " <<
+    // std::endl;
+
+    H5Dclose(_dset_id);
+    H5Sclose(_filespace);
+    H5Sclose(_memspace);
+    H5Pclose(_plist_id);
+}
+
+} // end namespace adios