Here is a minimal self contained example that produces the error.

On Thu, Sep 18, 2025 at 8:16 AM 'Wolfgang Bangerth' via deal.II User Group <
[email protected]> wrote:

> On 9/16/25 22:56, Jalil Khan wrote:
> > Here is how I call my data_out object to write HDF5,
> >
> > |DataOutBase::DataOutFilter
> > data_filter(DataOutBase::DataOutFilterFlags(true,true));
> > data_out.write_filtered_data(data_filter);
> > data_out.write_hdf5_parallel(data_filter, write_mesh_file,
> mesh_filename,
> > solution_filename, mpi_comm);|
>
> Jalil:
> Can you make this into a self-contained testcase that illustrates the
> problem?
> Best
>   W.
>
> --
> The deal.II project is located at http://www.dealii.org/
> For mailing list/forum options, see
> https://groups.google.com/d/forum/dealii?hl=en
> ---
> You received this message because you are subscribed to the Google Groups
> "deal.II User Group" group.
> To unsubscribe from this group and stop receiving emails from it, send an
> email to [email protected].
> To view this discussion visit
> https://groups.google.com/d/msgid/dealii/b5366d3f-6555-4683-9085-8ffd89c41218%40colostate.edu
> .
>


-- 
Regards,
JRK

-- 
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see 
https://groups.google.com/d/forum/dealii?hl=en
--- 
You received this message because you are subscribed to the Google Groups 
"deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to [email protected].
To view this discussion visit 
https://groups.google.com/d/msgid/dealii/CAALjy%3DgMG98fswUjoE4T7ffgm%3DRU_yhsXQp%3D9M8TM3HxkgxuRw%40mail.gmail.com.
#include <deal.II/base/utilities.h>
#include <deal.II/base/mpi.h>
#include <deal.II/distributed/tria.h>
#include <deal.II/grid/grid_generator.h>
#include <deal.II/fe/fe_q.h>
#include <deal.II/fe/fe_system.h>
#include <deal.II/dofs/dof_handler.h>
#include <deal.II/lac/la_parallel_vector.h>
#include <deal.II/numerics/data_out.h>
#include <deal.II/numerics/data_out_dof_data.h>
#include <deal.II/dofs/dof_tools.h>

#include <random>
#include <iostream>
#include <fstream>
#include <cstdlib>

using namespace dealii;

static const int dim = 3;
static const int degree = 4;

typedef LinearAlgebra::distributed::Vector<double> PVector;
typedef parallel::distributed::Triangulation<dim> PTriangulation;

//------------------------------------------------------------------------------
int main(int argc, char **argv)
{
  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);

  PTriangulation triangulation(MPI_COMM_WORLD);
  GridGenerator::hyper_cube(triangulation, 0.0, 1.0);
  triangulation.refine_global(7); // 2^7 x degree(4) approx 512
  DoFHandler<dim>  dof_handler(triangulation);
  const FESystem<dim> fe(FE_Q<dim>(degree), dim+2);
  dof_handler.clear();
  dof_handler.distribute_dofs(fe);
  const auto& locally_owned_dofs = dof_handler.locally_owned_dofs();
  IndexSet locally_relevant_dofs = DoFTools::extract_locally_relevant_dofs(dof_handler);
  PVector solution;
  solution.reinit(locally_owned_dofs, locally_relevant_dofs, MPI_COMM_WORLD);
  // Fill random values into solution
  std::mt19937 rng(42);
  std::uniform_real_distribution<double> dist(0.0, 1.0);
  for(auto i : solution.locally_owned_elements())
    solution[i] = dist(rng);
  solution.compress(VectorOperation::insert);
  solution.update_ghost_values();
  solution *= 0.1234;

  const bool write_mesh_file = true;
  std::string mesh_filename = "mesh.h5";
  std::string solution_filename = "vars.h5";

  DataOut<dim> data_out;
  data_out.add_data_vector(dof_handler, solution, "solution");
  data_out.build_patches(fe.degree);

  DataOutBase::DataOutFilter data_filter(DataOutBase::DataOutFilterFlags(true, true));
  data_out.write_filtered_data(data_filter);

  data_out.write_hdf5_parallel(data_filter,
                               write_mesh_file,
                               mesh_filename,
                               solution_filename,
                               MPI_COMM_WORLD);

  XDMFEntry new_xdmf_entry = data_out.create_xdmf_entry(data_filter,
                                                        mesh_filename,
                                                        solution_filename,
                                                        0.0,
                                                        MPI_COMM_WORLD);

  std::vector<XDMFEntry>        xdmf_entries;
  xdmf_entries.push_back(new_xdmf_entry);
  data_out.write_xdmf_file(xdmf_entries, "solution.xdmf", MPI_COMM_WORLD);

  return 0;
}

Reply via email to