From 1219b7f9d47b2cae64faaaca726c409d53532f46 Mon Sep 17 00:00:00 2001 From: Thomas Fischer <thomas.fischer@ufz.de> Date: Tue, 7 Jun 2022 13:13:52 +0200 Subject: [PATCH] [MeL/IO/VtkIO] Serial output function if mpi size is one. --- MeshLib/IO/VtkIO/PVDFile.cpp | 19 ++++++++++++++++--- MeshLib/IO/VtkIO/VtuInterface.cpp | 8 ++++++-- 2 files changed, 22 insertions(+), 5 deletions(-) diff --git a/MeshLib/IO/VtkIO/PVDFile.cpp b/MeshLib/IO/VtkIO/PVDFile.cpp index 2150b5caff9..e8fa6b26faa 100644 --- a/MeshLib/IO/VtkIO/PVDFile.cpp +++ b/MeshLib/IO/VtkIO/PVDFile.cpp @@ -14,6 +14,10 @@ #include <iomanip> #include <limits> +#ifdef USE_PETSC +#include <mpi.h> +#endif + #include "BaseLib/Error.h" #include "MeshLib/IO/VtkIO/VtuInterface.h" @@ -24,10 +28,19 @@ namespace IO void PVDFile::addVTUFile(const std::string& vtu_fname, double timestep) { #ifdef USE_PETSC - auto const vtu_file_name = - getVtuFileNameForPetscOutputWithoutExtension(vtu_fname); + int mpi_size; + MPI_Comm_size(MPI_COMM_WORLD, &mpi_size); + if (mpi_size == 1) + { + _datasets.emplace_back(timestep, vtu_fname); + } + else + { + auto const vtu_file_name = + getVtuFileNameForPetscOutputWithoutExtension(vtu_fname); - _datasets.emplace_back(timestep, vtu_file_name + ".pvtu"); + _datasets.emplace_back(timestep, vtu_file_name + ".pvtu"); + } #else _datasets.emplace_back(timestep, vtu_fname); #endif diff --git a/MeshLib/IO/VtkIO/VtuInterface.cpp b/MeshLib/IO/VtkIO/VtuInterface.cpp index 9c8615ee196..56a3b98c8b3 100644 --- a/MeshLib/IO/VtkIO/VtuInterface.cpp +++ b/MeshLib/IO/VtkIO/VtuInterface.cpp @@ -142,12 +142,16 @@ bool VtuInterface::writeToFile(std::filesystem::path const& file_path) MPI_Initialized(&mpi_init); if (mpi_init == 1) { + int mpi_size; + MPI_Comm_size(MPI_COMM_WORLD, &mpi_size); + if (mpi_size == 1) + { + return writeVTU<vtkXMLUnstructuredGridWriter>(file_path.string()); + } auto const vtu_file_name = getVtuFileNameForPetscOutputWithoutExtension(file_path.string()); int rank; MPI_Comm_rank(MPI_COMM_WORLD, &rank); - int mpi_size; - MPI_Comm_size(MPI_COMM_WORLD, &mpi_size); return writeVTU<vtkXMLPUnstructuredGridWriter>(vtu_file_name + ".pvtu", mpi_size, rank); } -- GitLab