Commit 96babbb6 authored by Jakub Klinkovský's avatar Jakub Klinkovský
Browse files

added options for tetrahedral remeshing

closes #1
parent 2104a4cf
Loading
Loading
Loading
Loading
Loading
+71 −4
Original line number Diff line number Diff line
@@ -13,6 +13,7 @@
#include <CGAL/Polyhedral_mesh_domain_with_features_3.h>
#include <CGAL/make_mesh_3.h>
#include <CGAL/refine_mesh_3.h>
#include <CGAL/tetrahedral_remeshing.h>
#include <CGAL/IO/output_to_vtu.h>

#include "lyra/lyra.hpp"
@@ -133,6 +134,21 @@ struct Parameters
	// degrees) of mesh cells.
	double exude_sliver_bound = 0;

// Multi-Material Isotropic Tetrahedral Remeshing
// - https://doc.cgal.org/latest/Tetrahedral_remeshing/index.html
// - https://doc.cgal.org/latest/Tetrahedral_remeshing/group__PkgTetrahedralRemeshingRef.html
	// remesh_iterations: The number of iterations for the remeshing algorithm.
	// Value of 0 disables remeshing.
	std::size_t remesh_iterations = 0;
	// remesh_target_edge_length: This parameter provides a uniform mesh
	// density target for the remeshing algorithm.
	double remesh_target_edge_length = 0;
	// remesh_boundaries: If false, none of the input volume boundaries can be
	// modified. Otherwise, the topology is preserved, but atomic operations
	// can be performed on the surfaces, and along feature polylines, such that
	// boundaries are remeshed.
	bool remesh_boundaries = true;

	void set_cli_options(lyra::cli & cli)
	{
		// pipeline options
@@ -284,6 +300,24 @@ struct Parameters
				["--exude-sliver-bound"]
				("A targeted lower bound on dihedral angles (in degrees) of mesh "
				 "cells." + fmt_default(exude_sliver_bound));

		// remeshing parameters
		cli |= lyra::opt(remesh_iterations, "N")
				["--remesh-iterations"]
				("The number of iterations for the Multi-Material Isotropic "
				 "Tetrahedral Remeshing algorithm. Value of 0 disables "
				 "remeshing." + fmt_default(remesh_iterations));
		cli |= lyra::opt(remesh_target_edge_length, "SCALAR")
				["--remesh-target-edge-length"]
				("This parameter provides a uniform mesh density target for the "
				 "remeshing algorithm." + fmt_default(remesh_target_edge_length));
		cli |= lyra::opt(remesh_boundaries, "BOOL")
				["--remesh-boundaries"]
				("If false, none of the input volume boundaries can be "
				 "modified. Otherwise, the topology is preserved, but atomic "
				 "operations can be performed on the surfaces, and along feature "
				 "polylines, such that boundaries are remeshed."
				 + fmt_default(exude_time_limit));
	}

	void print_options(std::ostream & str)
@@ -326,6 +360,11 @@ struct Parameters
			str << "  - Time limit: " << exude_time_limit << " seconds\n";
			str << "  - Sliver bound: " << exude_sliver_bound << \n";
		}
		str << "- Number of remeshing iterations: " << remesh_iterations << "\n";
		if (remesh_iterations > 0) {
			str << "  - Remeshing target edge length: " << remesh_target_edge_length << "\n";
			str << "  - Remesh boundaries: " << remesh_boundaries << "\n";
		}
		str.flush();
	}
};
@@ -374,6 +413,8 @@ void polyDomain(const std::filesystem::path & input_file,
	if (args.cell_size < 0) throw std::invalid_argument("cell_size cannot be negative");
	if (args.detect_features && args.edge_size == 0)
		throw std::invalid_argument("edge_size must be positive when detection of sharp features is enabled");
	if (args.remesh_iterations > 0 && args.remesh_target_edge_length <= 0)
		throw std::invalid_argument("remesh_iterations is positive, but remesh_target_edge_length is not positive");

	Polyhedron casePoly;
	readPolyhedron(casePoly, input_file);
@@ -436,11 +477,37 @@ void polyDomain(const std::filesystem::path & input_file,
			);
	}

	// Remeshing
	if (args.remesh_iterations > 0) {
		std::cout << "Remeshing..." << std::endl;
		// Extract triangulation from c3t3 and discard the c3t3 object
		using Triangulation_3 = CGAL::Triangulation_3<typename Tr::Geom_traits, typename Tr::Triangulation_data_structure>;
		Triangulation_3 t3 = CGAL::convert_to_triangulation_3(std::move(c3t3));

		// Remesh
		CGAL::tetrahedral_isotropic_remeshing(t3,
			args.remesh_target_edge_length,
			CGAL::parameters::number_of_iterations(args.remesh_iterations)
			.remesh_boundaries(args.remesh_boundaries)
		);

		// Create another c3t3 for output (but of different type)
		CGAL::Mesh_complex_3_in_triangulation_3<Triangulation_3, int, int> c3t3;
		c3t3.triangulation() = std::move(t3);
		c3t3.rescan_after_load_of_triangulation();

		// Output
		std::ofstream out(output_file);
		output_to_vtu(out, c3t3);
		out.close();
	}
	else {
		// Output
		std::ofstream out(output_file);
		output_to_vtu(out, c3t3);
		out.close();
	}
}

int main(int argc, char * argv[])
{