@@ -643,4 +643,138 @@ std::vector<U> read_data(std::string point_or_cell, std::string filename,
643643
644644 return values;
645645}
646+
647+ // / @brief Read a CG1 function from a VTKHDF format file.
648+ // /
649+ // / @tparam U Scalar type of mesh
650+ // / @param filename Name of the file to read from.
651+ // / @param mesh Mesh previously read from the same file.
652+ // / @param timestep The time step to read for time-dependent data.
653+ // / @return The function read from file.
654+ // /
655+ // / @note This only supports meshes with a single cell type as of now.
656+ template <std::floating_point U>
657+ fem::Function<U> read_CG1_function (std::string filename,
658+ std::shared_ptr<mesh::Mesh<U>> mesh,
659+ std::int32_t timestep = 0 )
660+ {
661+ auto element = basix::create_element<U>(
662+ basix::element::family::P,
663+ mesh::cell_type_to_basix_type (mesh->topology ()->cell_type ()), 1 ,
664+ basix::element::lagrange_variant::unset,
665+ basix::element::dpc_variant::unset, false );
666+
667+ hid_t h5file = hdf5::open_file (mesh->comm (), filename, " r" , true );
668+ std::string dataset_name = " /VTKHDF/PointData/u" ;
669+
670+ std::vector<std::int64_t > shape
671+ = hdf5::get_dataset_shape (h5file, dataset_name);
672+ hdf5::close_file (h5file);
673+
674+ const std::size_t block_size = shape[1 ];
675+
676+ auto P1
677+ = std::make_shared<fem::FunctionSpace<U>>(fem::create_functionspace<U>(
678+ mesh, std::make_shared<fem::FiniteElement<U>>(
679+ element, std::vector<std::size_t >{block_size})));
680+
681+ fem::Function<U> u_in (P1);
682+
683+ std::shared_ptr<const common::IndexMap> index_map (
684+ mesh->geometry ().index_map ());
685+
686+ std::int64_t range0 = index_map->local_range ()[0 ];
687+
688+ int npoints = index_map->size_local ();
689+
690+ std::array<std::int64_t , 2 > range{range0, range0 + npoints};
691+
692+ const auto values
693+ = io::VTKHDF::read_data (" Point" , filename, *mesh, range, timestep);
694+
695+ // Parallel distribution. For vector functions we distribute each component
696+ // separately.
697+ std::vector<std::vector<double >> values_s (block_size);
698+ for (auto & v : values_s)
699+ {
700+ v.reserve (values.size () / block_size);
701+ }
702+ for (std::size_t i = 0 ; i < values.size () / block_size; ++i)
703+ {
704+ for (int j = 0 ; j < block_size; ++j)
705+ {
706+ values_s[j].push_back (values[block_size * i + j]);
707+ }
708+ }
709+
710+ std::vector<std::int64_t > entities (range[1 ] - range[0 ]);
711+ std::iota (entities.begin (), entities.end (), range[0 ]);
712+
713+ MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
714+ const std::int64_t ,
715+ MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t , 2 >>
716+ entities_span (entities.data (),
717+ std::array<std::size_t , 2 >{entities.size (), 1 });
718+
719+ std::vector<std::pair<std::vector<std::int32_t >, std::vector<double >>>
720+ entities_values;
721+
722+ for (std::size_t i = 0 ; i < values_s.size (); ++i)
723+ {
724+ entities_values.push_back (dolfinx::io::distribute_entity_data<double >(
725+ *mesh->topology (), mesh->geometry ().input_global_indices (),
726+ mesh->geometry ().index_map ()->size_global (),
727+ mesh->geometry ().cmap ().create_dof_layout (), mesh->geometry ().dofmap (),
728+ mesh::cell_dim (mesh::CellType::point), entities_span, values_s[i]));
729+ }
730+
731+ auto num_vertices_per_cell
732+ = dolfinx::mesh::num_cell_vertices (mesh->topology ()->cell_type ());
733+ std::vector<std::int32_t > local_vertex_map (num_vertices_per_cell);
734+
735+ for (int i = 0 ; i < num_vertices_per_cell; ++i)
736+ {
737+ const auto v_to_d
738+ = u_in.function_space ()->dofmap ()->element_dof_layout ().entity_dofs (0 ,
739+ i);
740+ assert (v_to_d.size () == 1 );
741+ local_vertex_map[i] = v_to_d.front ();
742+ }
743+
744+ const auto tdim = mesh->topology ()->dim ();
745+ const auto c_to_v = mesh->topology ()->connectivity (tdim, 0 );
746+ std::vector<std::int32_t > vertex_to_dofmap (
747+ mesh->topology ()->index_map (0 )->size_local ()
748+ + mesh->topology ()->index_map (0 )->num_ghosts ());
749+
750+ for (int i = 0 ; i < mesh->topology ()->index_map (tdim)->size_local (); ++i)
751+ {
752+ const auto local_vertices = c_to_v->links (i);
753+ const auto local_dofs = u_in.function_space ()->dofmap ()->cell_dofs (i);
754+ for (int j = 0 ; j < num_vertices_per_cell; ++j)
755+ {
756+ vertex_to_dofmap[local_vertices[j]] = local_dofs[local_vertex_map[j]];
757+ }
758+ }
759+
760+ /*
761+ * After the data is read and distributed, we need to place the
762+ * retrieved values in the correct position in the function's array,
763+ * reading values and positions from `entities_values`.
764+ */
765+ for (std::size_t i = 0 ; i < entities_values[0 ].first .size (); ++i)
766+ {
767+ for (std::size_t j = 0 ; j < block_size; ++j)
768+ {
769+ u_in.x ()
770+ ->array ()[block_size * vertex_to_dofmap[entities_values[0 ].first [i]]
771+ + j]
772+ = entities_values[j].second [i];
773+ }
774+ }
775+
776+ u_in.x ()->scatter_fwd ();
777+
778+ return u_in;
779+ }
646780} // namespace dolfinx::io::VTKHDF
0 commit comments