cgal icon indicating copy to clipboard operation
cgal copied to clipboard

Meshing 3D image but see only one object

Open CheLamVien opened this issue 1 month ago • 10 comments

Hello, I am using CGAL 6.1.0 (Mesh_3) to mesh a labeled 3D image. My input image contains multiple labeled objects, and I expect each connected component (label > 0) to be meshed.

CGAL reports during initialization:

Searching for connected components...
  28 components were found.
Construct initial points...

So CGAL detects 28 components in the labelled image. I have already reduce the meshsize to an extremely small value (sizing_scale=0.0001).

However, in the resulting mesh output (.mesh), only one object appears. I have attached the input image (download here) and the dumped C3t3 (CGAL::dump_c3t3(c3t3, "out")) for inspection. My code:

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Mesh_triangulation_3.h>
#include <CGAL/Mesh_complex_3_in_triangulation_3.h>
#include <CGAL/Mesh_criteria_3.h>
#include <CGAL/Mesh_3/generate_label_weights.h>
#include <CGAL/Mesh_3/Detect_features_in_image.h>
#include <CGAL/Labeled_mesh_domain_3.h>
#include <CGAL/make_mesh_3.h>
#include <CGAL/Image_3.h>
#include <CGAL/IO/output_to_vtu.h>
#include <CGAL/Mesh_3/Construct_initial_points_labeled_image.h>
#include <CGAL/facets_in_complex_3_to_triangle_mesh.h>
#include <oneapi/tbb/global_control.h>
#include <filesystem>
#include <iostream>
#include <fstream>
#include <stdexcept>
#include <CGAL/Random.h>

using K = CGAL::Exact_predicates_inexact_constructions_kernel;
using Image_domain = CGAL::Labeled_mesh_domain_3<K>;
using Mesh_domain = CGAL::Mesh_domain_with_polyline_features_3<Image_domain>;

#ifdef CGAL_CONCURRENT_MESH_3
using Concurrency_tag = CGAL::Parallel_tag;
#else
using Concurrency_tag = CGAL::Sequential_tag;
#endif

using Tr = CGAL::Mesh_triangulation_3<Mesh_domain, CGAL::Default, Concurrency_tag>::type;
using C3t3 = CGAL::Mesh_complex_3_in_triangulation_3<Tr>;
using Mesh_criteria = CGAL::Mesh_criteria_3<Tr>;

namespace fs = std::filesystem;
using namespace CGAL::parameters;

int main(int argc, char* argv[])
{
  if (argc < 3) {
    throw std::invalid_argument("Usage: <program> <input_file.inr.gz> <output_folder> [sizing_scale] [edge_size]");
  }

  const std::string fname = argv[1];
  const std::string output_folder = argv[2];
  double sizing_scale = (argc > 3) ? std::atof(argv[3]) : 0.1;
  double edge_size = (argc > 4) ? std::atof(argv[4]) : 1000.0;

  CGAL::Image_3 image;
  if(!image.read(fname)){
    std::cerr << "Error: Cannot read file " << fname << std::endl;
    return EXIT_FAILURE;
  }

  const float sigma = (std::max)(image.vx(), (std::max)(image.vy(), image.vz()));
  CGAL::Image_3 img_weights = CGAL::Mesh_3::generate_label_weights(image, sigma);

  Mesh_domain domain = Mesh_domain::create_labeled_image_mesh_domain(
    image,
    weights = std::ref(img_weights),
    relative_error_bound = 1e-6,
    features_detector = CGAL::Mesh_3::Detect_features_in_image()
  );

  CGAL::Bbox_3 bbox = domain.bbox();
  double diag = CGAL::sqrt(
      CGAL::square(bbox.xmax() - bbox.xmin()) +
      CGAL::square(bbox.ymax() - bbox.ymin()) +
      CGAL::square(bbox.zmax() - bbox.zmin()));
  double sizing_default = diag * sizing_scale;

  Mesh_criteria criteria(
    edge_size = sizing_default,
    edge_distance = sizing_default / 10,
    facet_angle = 25,
    facet_size = sizing_default,
    facet_distance = sizing_default / 10,
    cell_radius_edge_ratio = 0,
    cell_size = 0
  );

  CGAL::Construct_initial_points_labeled_image<C3t3, Mesh_domain>
      img_pts_generator(image, domain);

  C3t3 c3t3 = CGAL::make_mesh_3<C3t3>(
    domain,
    criteria,
    initial_points_generator(img_pts_generator)
  );

  c3t3.remove_isolated_vertices();

  // Write MEDIT
  std::string output_file = output_folder + "/" + fs::path(fname).stem().string() + ".mesh";
  std::ofstream medit_file(output_file);
  CGAL::IO::write_MEDIT(medit_file, c3t3,
      all_cells(false).all_vertices(false).show_patches(false));
  medit_file.close();

  // Debug dump
  CGAL::dump_c3t3(c3t3, "out");

  return 0;
}

Any clarification about how Mesh_3 handles multiple connected labelled components would be very appreciated. If I am using the wrong export parameters or mesh criteria, please let me know.

Environment

  • Operating system: Ubuntu
  • CGAL version: 6.1.0

CheLamVien avatar Nov 28 '25 17:11 CheLamVien

Reproduced on Linux (GCC 13.3.0).

I was able to compile and run the reproduction code (using the provided input.inr.gz). The process completes successfully, but the resulting mesh is extremely small, confirming that most components are missing:

  • Vertices: 27
  • Cells: 119

Given that the input is expected to contain 28 labeled components, this low element count confirms that the output mesh likely contains only a single component.

LOTUS9679 avatar Dec 03 '25 09:12 LOTUS9679

As documented here, you can use something like in this example.

sloriot avatar Dec 03 '25 10:12 sloriot

Hi @sloriot I have already tried it in the given script, that's why I got the input "28 conponents were found". I would also like to note that all the objects are almost the same size, but only one can be meshed. I also have an example where the objects are much larger, but there are still many objects which can not be meshed.

Cheers, Vien Che

CheLamVien avatar Dec 03 '25 11:12 CheLamVien

Strange, I gave it a try in the CGAL lab demo and I could get all connected components. I'll have another look tomorrow directly with the code.

sloriot avatar Dec 03 '25 15:12 sloriot

Hi @sloriot, did you have a chance to review this issue again? I'm looking forward to hearing from you.

CheLamVien avatar Dec 05 '25 14:12 CheLamVien

Indeed it does not work.

@lrineau @janetournois I simply updated the CGAL example to add read from image: https://gist.github.com/sloriot/b4812ffda6c0c8d9780b0d6b8409258e

Using the data shared by @CheLamVien here

It does print that it found 28 components but only one is present in the output.

sloriot avatar Dec 10 '25 08:12 sloriot

Hello @CheLamVien , you may have seen that Laurent has provided a fix in PR #9178 . Let us know if it works for you.

afabri avatar Dec 10 '25 15:12 afabri

Hallo @afabri , I tested the fix provided by @lrineau. It works really well now. Thank you so much for your time.

Best regards, Vien Che

CheLamVien avatar Dec 10 '25 17:12 CheLamVien

Out of curiosity, what is it that you are meshing, and what is the use case?

afabri avatar Dec 11 '25 10:12 afabri

Hi @afabri, I meshed the cells from the fluorescent image of the cartilage tissue. I would like to conduct an FEM simulation for electrical stimulation in cartilage regeneration.

CheLamVien avatar Dec 11 '25 11:12 CheLamVien