Loading [MathJax]/extensions/tex2jax.js
OGS
FindCellsForPoints.h
Go to the documentation of this file.
1
10#pragma once
11
12#include <vtkCellData.h>
13#include <vtkCellTreeLocator.h>
14#include <vtkCharArray.h>
15#include <vtkExtractCells.h>
16#if VTK_VERSION_STRIPPED < 940
17#include <vtkIdFilter.h>
18#else
19#include <vtkGenerateIds.h>
20#endif
21#include <vtkIdList.h>
22#include <vtkNew.h>
23#include <vtkUnstructuredGrid.h>
24
25#include <Eigen/Core>
26#include <array>
27#include <range/v3/view/iota.hpp>
28#include <range/v3/view/transform.hpp>
29#include <vector>
30
31#include "BaseLib/Logging.h"
32
33namespace ApplicationUtils
34{
44{
45public:
46 void initialize(vtkUnstructuredGrid* const bulk_mesh)
47 {
48#if VTK_VERSION_STRIPPED < 940
49 vtkNew<vtkIdFilter> id_filter;
50 id_filter->SetInputData(bulk_mesh);
51 id_filter->SetCellIdsArrayName("bulk_element_ids");
52 id_filter->SetPointIds(false);
53 id_filter->SetCellIds(true);
54 id_filter->Update();
55 bulk_mesh_ = id_filter->GetOutput();
56#else
57 vtkNew<vtkGenerateIds> generate_ids;
58 generate_ids->SetInputData(bulk_mesh);
59 generate_ids->SetCellIdsArrayName("bulk_element_ids");
60 generate_ids->PointIdsOff();
61 generate_ids->CellIdsOn();
62 generate_ids->Update();
64 static_cast<vtkUnstructuredGrid*>(generate_ids->GetOutput());
65#endif
66 locator_ = vtkSmartPointer<vtkCellTreeLocator>::New();
67 locator_->SetDataSet(bulk_mesh_);
68 locator_->BuildLocator();
69 }
70
71 std::vector<vtkIdType> findCellIds(Eigen::Vector3d const& coords,
72 double const tolerance)
73 {
74 using namespace ranges;
75 auto const cell_ids = findCellIdsImpl(coords, tolerance);
76 return views::iota(0, cell_ids->GetNumberOfIds()) |
77 views::transform([&](vtkIdType i)
78 { return cell_ids->GetId(i); }) |
79 to_vector;
80 }
81
82 vtkSmartPointer<vtkUnstructuredGrid> find(Eigen::Vector3d const& coords,
83 double const tolerance)
84 {
85 auto const cell_ids = findCellIdsImpl(coords, tolerance);
86
87 vtkNew<vtkExtractCells> extract_cells;
88 extract_cells->SetInputData(bulk_mesh_);
89 extract_cells->SetCellList(cell_ids);
90 extract_cells->Update();
91 auto const filtered_bulk_mesh = extract_cells->GetOutput();
92
93 return filtered_bulk_mesh;
94 }
95
96private:
97 vtkSmartPointer<vtkIdList> findCellIdsImpl(Eigen::Vector3d const& coords,
98 double const tolerance)
99 {
100 std::array<double, 6> bbox = {
101 coords[0] - tolerance, coords[0] + tolerance,
102 coords[1] - tolerance, coords[1] + tolerance,
103 coords[2] - tolerance, coords[2] + tolerance};
104
105 vtkNew<vtkIdList> cell_ids;
106 locator_->FindCellsWithinBounds(bbox.data(), cell_ids);
107
108 return cell_ids;
109 }
110
111 vtkSmartPointer<vtkDataSet> bulk_mesh_;
112 vtkSmartPointer<vtkCellTreeLocator> locator_;
113};
114
115} // namespace ApplicationUtils
std::vector< vtkIdType > findCellIds(Eigen::Vector3d const &coords, double const tolerance)
vtkSmartPointer< vtkDataSet > bulk_mesh_
vtkSmartPointer< vtkIdList > findCellIdsImpl(Eigen::Vector3d const &coords, double const tolerance)
vtkSmartPointer< vtkCellTreeLocator > locator_
vtkSmartPointer< vtkUnstructuredGrid > find(Eigen::Vector3d const &coords, double const tolerance)
void initialize(vtkUnstructuredGrid *const bulk_mesh)