SemiDiscreteOT 1.0
Semi-Discrete Optimal Transport Library
Loading...
Searching...
No Matches
MeshManager.h
Go to the documentation of this file.
1#ifndef MESH_MANAGER_H
2#define MESH_MANAGER_H
3
4#include <deal.II/grid/tria.h>
5#include <deal.II/distributed/tria.h>
6#include <deal.II/distributed/fully_distributed_tria.h>
7#include <deal.II/grid/grid_generator.h>
8#include <deal.II/grid/grid_out.h>
9#include <deal.II/grid/grid_in.h>
10#include <deal.II/grid/grid_tools.h>
11#include <deal.II/base/utilities.h>
12#include <deal.II/base/conditional_ostream.h>
13#include <deal.II/numerics/data_out.h>
14
15#include <filesystem>
16#include <vector>
17#include <string>
18#include <fstream>
19#include <stdexcept>
20
21using namespace dealii;
22
23template <int dim, int spacedim = dim>
25public:
30 MeshManager(const MPI_Comm& comm);
31
40 template<typename TriangulationType>
41 void generate_mesh(TriangulationType& tria,
42 const std::string& grid_generator_function,
43 const std::string& grid_generator_arguments,
44 const unsigned int n_refinements,
45 const bool use_tetrahedral_mesh);
46
52 void load_source_mesh(parallel::fullydistributed::Triangulation<dim, spacedim>& source_mesh,
53 const std::string& source_file = "");
54
59 void load_target_mesh(Triangulation<dim, spacedim>& target_mesh);
60
67 void load_mesh_at_level(parallel::fullydistributed::Triangulation<dim, spacedim>& source_mesh,
68 DoFHandler<dim, spacedim>& dof_handler_source,
69 const std::string& mesh_file);
70
76 void save_meshes(const parallel::fullydistributed::Triangulation<dim, spacedim>& source_mesh,
77 const Triangulation<dim, spacedim>& target_mesh);
78
88 template<typename TriangulationType>
89 bool write_mesh(const TriangulationType& mesh,
90 const std::string& filepath,
91 const std::vector<std::string>& formats,
92 const std::vector<double>* cell_data = nullptr,
93 const std::string& data_name = "cell_data");
94
100 static std::vector<std::string> get_mesh_hierarchy_files(
101 const std::string& dir = "output/data_multilevel/source_multilevel");
102
103private:
105 const unsigned int n_mpi_processes;
106 const unsigned int this_mpi_process;
107 ConditionalOStream pcout;
108
109 // Constants
110 const std::string mesh_directory = "output/data_mesh";
111};
112
113// Include template implementations
115
116#endif // MESH_MANAGER_H
static std::vector< std::string > get_mesh_hierarchy_files(const std::string &dir="output/data_multilevel/source_multilevel")
Get sorted list of mesh hierarchy files.
MPI_Comm mpi_communicator
const unsigned int n_mpi_processes
bool write_mesh(const TriangulationType &mesh, const std::string &filepath, const std::vector< std::string > &formats, const std::vector< double > *cell_data=nullptr, const std::string &data_name="cell_data")
Write mesh to file in specified formats with optional cell data.
void load_source_mesh(parallel::fullydistributed::Triangulation< dim, spacedim > &source_mesh, const std::string &source_file="")
Load source mesh from file into distributed triangulation.
void load_target_mesh(Triangulation< dim, spacedim > &target_mesh)
Load target mesh from file into serial triangulation.
void generate_mesh(TriangulationType &tria, const std::string &grid_generator_function, const std::string &grid_generator_arguments, const unsigned int n_refinements, const bool use_tetrahedral_mesh)
Generate mesh using deal.II grid generator.
const std::string mesh_directory
const unsigned int this_mpi_process
ConditionalOStream pcout
void save_meshes(const parallel::fullydistributed::Triangulation< dim, spacedim > &source_mesh, const Triangulation< dim, spacedim > &target_mesh)
Save source and target meshes to files.
void load_mesh_at_level(parallel::fullydistributed::Triangulation< dim, spacedim > &source_mesh, DoFHandler< dim, spacedim > &dof_handler_source, const std::string &mesh_file)
Load mesh at specific refinement level.