12#include <tclap/CmdLine.h> 
   15#include <boost/property_tree/ptree.hpp> 
   16#include <boost/property_tree/xml_parser.hpp> 
   21#include <range/v3/algorithm/copy.hpp> 
   22#include <range/v3/algorithm/transform.hpp> 
   23#include <range/v3/numeric/accumulate.hpp> 
   24#include <range/v3/range/conversion.hpp> 
   25#include <range/v3/view/filter.hpp> 
   27#include <unordered_set> 
   57    std::vector<std::unique_ptr<MeshLib::Mesh>> 
const& partitioned_meshes,
 
   60    std::vector<MeshEntityMapInfo> 
const& merged_node_map,
 
   61    std::vector<MeshEntityMapInfo> 
const& merged_element_map)
 
   80    if (pv_name == 
"OGS_VERSION" || pv_name == 
"IntegrationPointMetaData")
 
   83            merged_mesh, pv_name, item_type, pv_num_components);
 
   89    std::vector<MeshLib::PropertyVector<T>*> partition_property_vectors;
 
   90    partition_property_vectors.reserve(partitioned_meshes.size());
 
   91    for (
auto const& mesh : partitioned_meshes)
 
   93        partition_property_vectors.emplace_back(
 
   94            mesh->getProperties().getPropertyVector<T>(pv_name, item_type,
 
   98    auto createNewCellOrNodePropertyVector =
 
   99        [&](std::vector<MeshEntityMapInfo> 
const& mesh_entity_map)
 
  102            merged_mesh, pv_name, item_type, pv_num_components);
 
  103        std::size_t counter = 0;
 
  104        for (
auto const& entity_info : mesh_entity_map)
 
  106            auto const& partition_pv =
 
  107                partition_property_vectors[entity_info.partition_id];
 
  108            for (
int i_com = 0; i_com < pv_num_components; i_com++)
 
  110                (*new_pv)[counter * pv_num_components + i_com] =
 
  111                    (*partition_pv)[entity_info.original_id *
 
  121        createNewCellOrNodePropertyVector(merged_node_map);
 
  127        createNewCellOrNodePropertyVector(merged_element_map);
 
  133        std::vector<std::vector<std::size_t>> partition_element_offsets;
 
  134        partition_element_offsets.reserve(partitioned_meshes.size());
 
  135        for (
auto const& mesh : partitioned_meshes)
 
  137            partition_element_offsets.emplace_back(
 
  139                    mesh->getElements(), *pv, properties));
 
  143            merged_mesh, pv_name, item_type, pv_num_components);
 
  146        auto const ip_meta_data =
 
  150        auto number_of_integration_points = [&](
auto const* element)
 
  153                ip_meta_data, *element);
 
  156        std::size_t 
const counter = ranges::accumulate(
 
  158                ranges::views::transform(number_of_integration_points),
 
  160        new_pv->resize(counter * pv_num_components);
 
  162        auto const global_ip_offsets =
 
  166        std::size_t element_counter = 0;
 
  167        for (
auto const& element_info : merged_element_map)
 
  170                *(partition_property_vectors[element_info.partition_id]);
 
  172            auto const& offsets =
 
  173                partition_element_offsets[element_info.partition_id];
 
  175            int const begin_pos = offsets[element_info.original_id];
 
  176            int const end_pos = offsets[element_info.original_id + 1];
 
  178            std::copy(partition_pv.
begin() + begin_pos,
 
  179                      partition_pv.
begin() + end_pos,
 
  180                      new_pv->begin() + global_ip_offsets[element_counter]);
 
 
  191    std::ifstream ins(pvtu_file_name);
 
  195        OGS_FATAL(
"Could not open pvtu file {:s}.", pvtu_file_name);
 
  198    using boost::property_tree::ptree;
 
  202    auto root = pt.get_child(
"VTKFile");
 
  204    std::vector<std::string> vtu_file_names;
 
  208    for (ptree::value_type 
const& v : root.get_child(
"PUnstructuredGrid"))
 
  210        if (v.first == 
"Piece")
 
  215                std::filesystem::path(v.second.get(
"<xmlattr>.Source", 
""))
 
  221    if (vtu_file_names.empty())
 
  223        OGS_FATAL(
"PVTU file {:s} does not contain any vtu piece",
 
  227    return vtu_file_names;
 
 
  231std::tuple<std::vector<MeshLib::Node*>, std::vector<std::size_t>>
 
  234    std::vector<std::size_t> number_of_nodes_per_partition;
 
  235    ranges::transform(meshes, std::back_inserter(number_of_nodes_per_partition),
 
  237                      { 
return mesh->getNumberOfNodes(); });
 
  238    std::vector<std::size_t> 
const offsets =
 
  241    std::vector<MeshLib::Node*> all_nodes;
 
  242    all_nodes.reserve(offsets.back());
 
  243    for (
auto const& mesh : meshes)
 
  245        ranges::copy(mesh->getNodes(), std::back_inserter(all_nodes));
 
  247    return {all_nodes, offsets};
 
 
  250std::tuple<std::vector<MeshLib::Element*>, std::vector<MeshEntityMapInfo>>
 
  253    std::vector<MeshLib::Element*> regular_elements;
 
  255    std::size_t partition_counter = 0;
 
  256    std::vector<MeshEntityMapInfo> merged_element_map;
 
  257    for (
auto const& mesh : meshes)
 
  261        auto const* 
const ghost_id_vector =
 
  263        assert(ghost_id_vector);
 
  265        auto const& mesh_elements = mesh->getElements();
 
  267        auto const last_element_id_of_previous_partition =
 
  268            regular_elements.size();
 
  271        { 
return (*ghost_id_vector)[element->getID()] == 0; };
 
  273        auto regular_mesh_elements =
 
  274            mesh_elements | ranges::views::filter(is_regular_element);
 
  275        regular_elements.insert(regular_elements.end(),
 
  276                                regular_mesh_elements.begin(),
 
  277                                regular_mesh_elements.end());
 
  279        for (
auto element_id = last_element_id_of_previous_partition;
 
  280             element_id < regular_elements.size();
 
  283            merged_element_map.push_back(
 
  284                {partition_counter, regular_elements[element_id]->getID()});
 
  290    return {regular_elements, merged_element_map};
 
 
  297    std::vector<MeshLib::Node*> query_nodes;
 
  299        std::numeric_limits<double>::epsilon() * extent.squaredNorm();
 
  301    auto constexpr dir = std::numeric_limits<double>::max();
 
  302    Eigen::Vector3d 
const max = {std::nextafter(node[0] + eps, dir),
 
  303                                 std::nextafter(node[1] + eps, dir),
 
  304                                 std::nextafter(node[2] + eps, dir)};
 
  307    if (query_nodes.empty())
 
  310            "query_nodes for node [{}], ({}, {}, {}) of element " 
  311            "[{}] are empty, eps is {}",
 
  312            node.
getID(), node[0], node[1], node[2], element_id, eps);
 
  315        std::find_if(query_nodes.begin(), query_nodes.end(),
 
  316                     [&node, eps](
auto const* p)
 
  318                         return (p->asEigenVector3d() - node.asEigenVector3d())
 
  319                                    .squaredNorm() < eps;
 
  321    if (it == query_nodes.end())
 
  324            "did not find node: [{}] ({}, {}, {}) of element [{}] in " 
  326            node.
getID(), node[0], node[1], node[2], element_id);
 
 
  332    std::vector<MeshLib::Element*> 
const& regular_elements,
 
  334    Eigen::Vector3d 
const& extent)
 
  336    for (
auto& e : regular_elements)
 
  338        for (
unsigned i = 0; i < e->getNumberOfNodes(); i++)
 
  340            auto* 
const node = e->getNode(i);
 
  342            if (!oct_tree.
addPoint(node, node_ptr))
 
  345                    oct_tree, extent, *node, e->getID());
 
  346                e->setNode(i, node_ptr);
 
 
  352std::pair<std::vector<MeshLib::Node*>, std::vector<MeshEntityMapInfo>>
 
  354                std::vector<std::size_t> 
const& partition_offsets,
 
  357    std::vector<MeshLib::Node*> unique_merged_nodes;
 
  358    unique_merged_nodes.reserve(all_merged_nodes_tmp.size());
 
  360    std::vector<MeshEntityMapInfo> merged_node_map;
 
  361    merged_node_map.reserve(all_merged_nodes_tmp.size());
 
  363    for (std::size_t i = 0; i < partition_offsets.size() - 1; ++i)
 
  365        for (std::size_t pos = partition_offsets[i];
 
  366             pos < partition_offsets[i + 1];
 
  369            auto* node = all_merged_nodes_tmp[pos];
 
  371            if (oct_tree.
addPoint(node, node_ptr))
 
  373                unique_merged_nodes.push_back(node);
 
  374                merged_node_map.push_back({i, pos - partition_offsets[i]});
 
  378    return {unique_merged_nodes, merged_node_map};
 
 
  382    std::vector<std::unique_ptr<MeshLib::Mesh>> 
const& meshes)
 
  385    merged_element_timer.
start();
 
  390        "Collection of {} regular elements and computing element map took {} s",
 
  391        regular_elements.size(), merged_element_timer.
elapsed());
 
  395    collect_nodes_timer.
start();
 
  396    auto [all_merged_nodes_tmp, partition_offsets] =
 
  398    INFO(
"Collection of {} nodes and computing offsets took {} s",
 
  399         all_merged_nodes_tmp.size(), collect_nodes_timer.
elapsed());
 
  402    merged_nodes_timer.
start();
 
  403    GeoLib::AABB aabb(all_merged_nodes_tmp.begin(), all_merged_nodes_tmp.end());
 
  404    auto oct_tree = std::unique_ptr<GeoLib::OctTree<MeshLib::Node, 16>>(
 
  408    auto [unique_merged_nodes, merged_node_map] =
 
  410    INFO(
"Make nodes unique ({} unique nodes) / computing map took {} s",
 
  411         unique_merged_nodes.size(), merged_nodes_timer.
elapsed());
 
  414    reset_nodes_in_elements_timer.
start();
 
  417    INFO(
"Reset nodes in regular elements took {} s",
 
  418         reset_nodes_in_elements_timer.
elapsed());
 
  421    mesh_creation_timer.
start();
 
  424    auto merged_mesh = std::make_unique<MeshLib::Mesh>(
 
  425        "pvtu_merged_mesh", unique_merged_nodes, regular_elements,
 
  427    INFO(
"creation of merged mesh took {} s", mesh_creation_timer.
elapsed());
 
  429    auto const& properties = meshes[0]->getProperties();
 
  432    property_timer.
start();
 
  435        [&, &merged_node_map = merged_node_map,
 
  436         &merged_element_map = merged_element_map](
auto type,
 
  437                                                   auto const& property)
 
  440                *merged_mesh, meshes,
 
  443                properties, merged_node_map, merged_element_map);
 
  445    INFO(
"merge properties into merged mesh took {} s",
 
 
  454              std::string_view 
const property_name,
 
  459        INFO(
"Skipping property '{}' since it is not from type '{}'.",
 
  460             property_name, 
typeid(T).name());
 
  463    auto const* subdomain_pv =
 
  465    auto const number_of_components =
 
  466        subdomain_pv->getNumberOfGlobalComponents();
 
  469        original_mesh, std::string(property_name), mesh_item_type,
 
  470        number_of_components);
 
  472    for (std::size_t i = 0; i < global_ids.
size(); ++i)
 
  474        auto const& global_id = global_ids[i];
 
  475        for (std::remove_cvref_t<
decltype(number_of_components)>
 
  476                 component_number = 0;
 
  477             component_number < number_of_components;
 
  480            original_pv->
getComponent(global_id, component_number) =
 
  481                subdomain_pv->getComponent(i, component_number);
 
  484    INFO(
"Data array {} from data type '{}' transferred.", property_name,
 
 
  493    auto const& subdomain_properties = subdomain_mesh.
getProperties();
 
  494    if (!subdomain_properties.existsPropertyVector<std::size_t>(
 
  498            "The data array '{}' is required for the transfer of data from " 
  499            "subdomain mesh to the global mesh. But it doesn't exist in " 
  500            "subdomain mesh '{}'.",
 
  503    auto const* global_ids =
 
  504        subdomain_properties.getPropertyVector<std::size_t>(
 
  506    if (global_ids == 
nullptr)
 
  508        OGS_FATAL(
"The data array '{}' is not available but required.",
 
  512    auto const& property_names =
 
  513        subdomain_properties.getPropertyVectorNames(mesh_item_type);
 
  515    for (
auto const& property_name : property_names)
 
  524                "Skipping property '{}' since it is adjusted locally in " 
  525                "NodeWiseMeshPartitioner::renumberBulkNodeIdsProperty().",
 
  529        if (property_name == 
"vtkGhostType")
 
  532                "Skipping property 'vtkGhostType' since it is only required " 
  533                "for parallel execution.");
 
  537                             property_name, mesh_item_type) ||
 
  539                            property_name, mesh_item_type) ||
 
  541                               property_name, mesh_item_type) ||
 
  543                                    *global_ids, property_name,
 
  546                                  *global_ids, property_name, mesh_item_type) ||
 
  547            transfer<int>(subdomain_properties, original_mesh, *global_ids,
 
  548                          property_name, mesh_item_type) ||
 
  550                           property_name, mesh_item_type))
 
  552            INFO(
"Data array {} transferred from '{}' to'{}'", property_name,
 
 
  568int main(
int argc, 
char* argv[])
 
  571        "This tool merges VTU files of PVTU into one single VTU file. Apart " 
  572        "from the mesh data, all property data are merged as well" 
  573        "\n\nOpenGeoSys-6 software, version " +
 
  576            "Copyright (c) 2012-2025, OpenGeoSys Community " 
  577            "(http://www.opengeosys.org)",
 
  580    TCLAP::ValueArg<std::string> output_arg(
 
  581        "o", 
"output", 
"Output (.vtu). The output mesh file (format = *.vtu)",
 
  582        true, 
"", 
"OUTPUT_FILE");
 
  585    TCLAP::ValueArg<std::string> input_arg(
 
  586        "i", 
"input", 
"Input (.pvtu). The partitioned input mesh file", 
true,
 
  590    TCLAP::ValueArg<std::string> original_mesh_input_arg(
 
  591        "m", 
"original_mesh",
 
  592        "optional, the original unpartitioned input mesh (*.vtu)", 
false, 
"",
 
  594    cmd.add(original_mesh_input_arg);
 
  598    cmd.add(log_level_arg);
 
  599    cmd.parse(argc, argv);
 
  606        OGS_FATAL(
"The extension of input file name {:s} is not \"pvtu\"",
 
  607                  input_arg.getValue());
 
  611        OGS_FATAL(
"The extension of output file name {:s} is not \"vtu\"",
 
  612                  output_arg.getValue());
 
  617    std::vector<std::unique_ptr<MeshLib::Mesh>> meshes;
 
  618    meshes.reserve(vtu_file_names.size());
 
  622    for (
auto const& file_name : vtu_file_names)
 
  624        auto mesh = std::unique_ptr<MeshLib::Mesh>(
 
  632                "Property vector vtkGhostType does not exist in mesh {:s}.",
 
  636        meshes.emplace_back(std::move(mesh));
 
  638    INFO(
"Reading meshes took {} s", io_timer.
elapsed());
 
  640    std::unique_ptr<MeshLib::Mesh> merged_mesh;
 
  641    if (original_mesh_input_arg.getValue().empty())
 
  651        merged_mesh = std::unique_ptr<MeshLib::Mesh>(
 
  653                original_mesh_input_arg.getValue()));
 
  654        INFO(
"Reading original unpartitioned mesh took {} s",
 
  656        for (
auto const& subdomain_mesh : meshes)
 
  659                *(subdomain_mesh.get()), *merged_mesh);
 
  666    writing_timer.
start();
 
  667    auto const result = writer.
writeToFile(output_arg.getValue());
 
  670        ERR(
"Could not write mesh to '{:s}'.", output_arg.getValue());
 
  673    INFO(
"writing mesh took {} s", writing_timer.
elapsed());
 
  680    merged_mesh->shallowClean();
 
 
Definition of the AABB class.
 
Definition of the Element class.
 
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
 
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
 
Definition of the class Properties that implements a container of properties.
 
Definition of the Mesh class.
 
Definition of the Node class.
 
Implementation of the OctTree class.
 
bool transfer(MeshLib::Properties const &subdomain_properties, MeshLib::Mesh &original_mesh, MeshLib::PropertyVector< std::size_t > const &global_ids, std::string_view const property_name, MeshLib::MeshItemType const mesh_item_type)
 
int main(int argc, char *argv[])
 
std::unique_ptr< MeshLib::Mesh > mergeSubdomainMeshes(std::vector< std::unique_ptr< MeshLib::Mesh > > const &meshes)
 
std::pair< std::vector< MeshLib::Node * >, std::vector< MeshEntityMapInfo > > makeNodesUnique(std::vector< MeshLib::Node * > const &all_merged_nodes_tmp, std::vector< std::size_t > const &partition_offsets, GeoLib::OctTree< MeshLib::Node, 16 > &oct_tree)
 
MeshLib::Node * getExistingNodeFromOctTree(GeoLib::OctTree< MeshLib::Node, 16 > &oct_tree, Eigen::Vector3d const &extent, MeshLib::Node const &node, std::size_t const element_id)
 
void transferPropertiesFromPartitionedMeshToUnpartitionedMesh(MeshLib::Mesh const &subdomain_mesh, MeshLib::Mesh &original_mesh, MeshLib::MeshItemType const mesh_item_type)
 
void resetNodesInRegularElements(std::vector< MeshLib::Element * > const ®ular_elements, GeoLib::OctTree< MeshLib::Node, 16 > &oct_tree, Eigen::Vector3d const &extent)
 
std::tuple< std::vector< MeshLib::Element * >, std::vector< MeshEntityMapInfo > > getRegularElements(std::vector< std::unique_ptr< MeshLib::Mesh > > const &meshes)
 
std::vector< std::string > readVtuFileNames(std::string const &pvtu_file_name)
 
bool createPropertyVector(MeshLib::Mesh &merged_mesh, std::vector< std::unique_ptr< MeshLib::Mesh > > const &partitioned_meshes, MeshLib::PropertyVector< T > const *const pv, MeshLib::Properties const &properties, std::vector< MeshEntityMapInfo > const &merged_node_map, std::vector< MeshEntityMapInfo > const &merged_element_map)
 
std::tuple< std::vector< MeshLib::Node * >, std::vector< std::size_t > > getMergedNodesVector(std::vector< std::unique_ptr< MeshLib::Mesh > > const &meshes)
 
Definition of the Point3d class.
 
void applyToPropertyVectors(Properties const &properties, Function f)
 
Definition of the RunTime class.
 
Implementation of the VtuInterface class.
 
double elapsed() const
Get the elapsed time in seconds.
 
void start()
Start the timer.
 
Class AABB is an axis aligned bounding box around a given set of geometric points of (template) type ...
 
Eigen::Vector3d const & getMaxPoint() const
 
Eigen::Vector3d const & getMinPoint() const
 
bool addPoint(POINT *pnt, POINT *&ret_pnt)
 
void getPointsInRange(T const &min, T const &max, std::vector< POINT * > &pnts) const
 
static OctTree< POINT, MAX_POINTS > * createOctTree(Eigen::Vector3d ll, Eigen::Vector3d ur, double eps=std::numeric_limits< double >::epsilon())
 
std::size_t getID() const
 
Eigen::Vector3d const & asEigenVector3d() const
 
Reads and writes VtkXMLUnstructuredGrid-files (vtu) to and from OGS data structures....
 
static MeshLib::Mesh * readVTUFile(std::string const &file_name, bool const compute_element_neighbors=false)
 
bool writeToFile(std::filesystem::path const &file_path)
 
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
 
Properties & getProperties()
 
const std::string getName() const
Get name of the mesh.
 
Property manager on mesh items. Class Properties manages scalar, vector or matrix properties....
 
bool existsPropertyVector(std::string_view name) const
 
PropertyVector< T > const * getPropertyVector(std::string_view name) const
 
MeshItemType getMeshItemType() const
 
int getNumberOfGlobalComponents() const
 
std::string const & getPropertyName() const
 
constexpr std::size_t size() const
 
PROP_VAL_TYPE & getComponent(std::size_t tuple_index, int component)
Returns the value for the given component stored in the given tuple.
 
constexpr PROP_VAL_TYPE * begin()
 
TCLAP::ValueArg< std::string > makeLogLevelArg()
 
void initOGSLogger(std::string const &log_level)
 
std::string getFileExtension(const std::string &path)
 
std::string extractPath(std::string const &pathname)
 
std::string joinPaths(std::string const &pathA, std::string const &pathB)
 
std::vector< ranges::range_value_t< R > > sizesToOffsets(R const &sizes)
 
GITINFOLIB_EXPORT const std::string ogs_version
 
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)
 
std::optional< IntegrationPointMetaData > getIntegrationPointMetaData(MeshLib::Properties const &properties)
 
constexpr std::string_view getBulkIDString(MeshItemType mesh_item_type)
 
constexpr std::string_view globalIDString(MeshLib::MeshItemType const mesh_item_type)
 
IntegrationPointMetaDataSingleField getIntegrationPointMetaDataSingleField(std::optional< IntegrationPointMetaData > const &ip_meta_data, std::string const &field_name)