Chaste  Build::
ImageToMesh.cpp
1 /*
2 
3 Copyright (c) 2005-2016, University of Oxford.
4  All rights reserved.
5 
6  University of Oxford means the Chancellor, Masters and Scholars of the
7  University of Oxford, having an administrative office at Wellington
8  Square, Oxford OX1 2JD, UK.
9 
10  This file is part of Chaste.
11 
12  Redistribution and use in source and binary forms, with or without
13  modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24  AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25  IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26  ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27  LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28  CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29  GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30  HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31  LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32  OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34  */
35 
36 #include "Exception.hpp"
37 #include "ImageToMesh.hpp"
38 #include <vtkThreshold.h>
39 #include <vtkGeometryFilter.h>
40 #include <vtkMarchingCubes.h>
41 #include <vtkTriangleFilter.h>
42 #include <vtkPoints.h>
43 #include <vtkPolyData.h>
44 #include <vtkProbeFilter.h>
45 #include "ImageToSurface.hpp"
46 #include "SurfaceCleaner.hpp"
47 #include "BoundaryExtractor.hpp"
48 #include "UblasIncludes.hpp"
49 #include "GeometryWriter.hpp"
50 #include "DiscreteContinuumMeshGenerator.hpp"
51 #include "BaseUnits.hpp"
52 
53 template<unsigned DIM>
55  : mpImage(vtkSmartPointer<vtkImageData>::New()),
56  mElementSize(0.0 * unit::metres * unit::metres* unit::metres),
57  mMesh(),
58  mpDomain()
59 {
60 
61 }
62 
63 template<unsigned DIM>
64 boost::shared_ptr<ImageToMesh<DIM> > ImageToMesh<DIM>::Create()
65 {
66  MAKE_PTR(ImageToMesh, pSelf);
67  return pSelf;
68 }
69 
70 template<unsigned DIM>
72 {
73 
74 }
75 
76 template<unsigned DIM>
77 boost::shared_ptr<DiscreteContinuumMesh<DIM, DIM> > ImageToMesh<DIM>::GetMesh()
78 {
79  if(mMesh)
80  {
81  return(mMesh);
82  }
83  else
84  {
85  EXCEPTION("No mesh set. Did you run 'Update()' ?");
86  }
87 }
88 
89 template<unsigned DIM>
90 void ImageToMesh<DIM>::SetTissueDomain(boost::shared_ptr<Part<DIM> > pTissueDomain)
91 {
92  mpDomain = pTissueDomain;
93 }
94 
95 template<unsigned DIM>
96 void ImageToMesh<DIM>::SetElementSize(units::quantity<unit::volume> elementSize)
97 {
98  mElementSize = elementSize;
99 }
100 
101 template<unsigned DIM>
102 void ImageToMesh<DIM>::SetInput(vtkSmartPointer<vtkImageData> pImage)
103 {
104  mpImage = pImage;
105 }
106 
107 template<unsigned DIM>
108 void ImageToMesh<DIM>::SetInputRaw(vtkImageData* pImage)
109 {
110  mpImage.TakeReference(pImage);
111 // pImage->Delete();
112 }
113 
114 template<unsigned DIM>
116 {
117  if(!mpImage)
118  {
119  EXCEPTION("No input set.");
120  }
121 
122  ImageToSurface image_to_surface;
123  image_to_surface.SetInput(mpImage);
124  image_to_surface.SetThreshold(1.0, false);
125  image_to_surface.Update();
126 
127  BoundaryExtractor extractor;
128  extractor.SetInput(image_to_surface.GetOutput());
129  extractor.SetDoSmoothing(true);
130  extractor.SetSmoothingLength(10.0);
131  extractor.Update();
132 
133  if(DIM==2)
134  {
135  DiscreteContinuumMeshGenerator<DIM, DIM> temp_mesh_generator;
136  temp_mesh_generator.SetDomain(extractor.GetOutput());
137  if(mpDomain)
138  {
139  temp_mesh_generator.SetDomain(mpDomain);
140  }
141  temp_mesh_generator.Update();
142  boost::shared_ptr<DiscreteContinuumMesh<DIM, DIM> > p_temp_mesh = temp_mesh_generator.GetMesh();
143 
144  std::vector<DimensionalChastePoint<DIM> > mesh_points = p_temp_mesh->GetNodeLocationsAsPoints();
145  std::vector<std::vector<unsigned> > mesh_connectivity = p_temp_mesh->GetConnectivity();
146 
147  units::quantity<unit::length> reference_length = p_temp_mesh->GetReferenceLengthScale();
148  vtkSmartPointer<vtkPoints> p_vtk_points = vtkSmartPointer<vtkPoints>::New();
149  for(unsigned idx = 0; idx<mesh_connectivity.size(); idx++)
150  {
151  c_vector<double, DIM> location0 = mesh_points[mesh_connectivity[idx][0]].GetLocation(reference_length);
152  c_vector<double, DIM> location1 = mesh_points[mesh_connectivity[idx][1]].GetLocation(reference_length);
153  c_vector<double, DIM> location2 = mesh_points[mesh_connectivity[idx][2]].GetLocation(reference_length);
154  double centx = (location0[0] + location1[0] + location2[0])/3.0;
155  double centy = (location0[1] + location1[1] + location2[1])/3.0;
156  p_vtk_points->InsertNextPoint(centx, centy, 0.0);
157  }
158 
159  // Get the value of the vessel from the image to determine if the point is inside or outside a vessel
160  vtkSmartPointer<vtkPolyData> p_sampling_data = vtkSmartPointer<vtkPolyData>::New();
161  p_sampling_data->SetPoints(p_vtk_points);
162 
163  vtkSmartPointer<vtkProbeFilter> p_image_probe = vtkSmartPointer<vtkProbeFilter>::New();
164  #if VTK_MAJOR_VERSION <= 5
165  p_image_probe->SetInput(p_sampling_data);
166  p_image_probe->SetSource(mpImage);
167  #else
168  p_image_probe->SetInputData(p_sampling_data);
169  p_image_probe->SetSourceData(mpImage);
170  #endif
171  p_image_probe->Update();
172 
173  std::vector<DimensionalChastePoint<DIM> > holes;
174  std::vector<DimensionalChastePoint<DIM> > regions;
175  for(unsigned idx=0; idx<p_image_probe->GetOutput()->GetPointData()->GetScalars()->GetNumberOfTuples(); idx++)
176  {
177  if(p_image_probe->GetOutput()->GetPointData()->GetScalars()->GetTuple1(idx)>5)
178  {
179  if(mpDomain)
180  {
181  c_vector<double, DIM> loc;
182  loc[0] = p_vtk_points->GetPoint(idx)[0];
183  loc[1] = p_vtk_points->GetPoint(idx)[1];
184  if(mpDomain->GetPolygons()[0]->ContainsPoint(DimensionalChastePoint<DIM>(loc, reference_length)))
185  {
186  regions.push_back(DimensionalChastePoint<DIM>(loc, reference_length));
187  }
188  else
189  {
190  holes.push_back(DimensionalChastePoint<DIM>(loc, reference_length));
191  }
192  }
193  else
194  {
195  c_vector<double, DIM> loc;
196  loc[0] = p_vtk_points->GetPoint(idx)[0];
197  loc[1] = p_vtk_points->GetPoint(idx)[1];
198  holes.push_back(DimensionalChastePoint<DIM>(loc, reference_length));
199  }
200  }
201  }
202  DiscreteContinuumMeshGenerator<DIM, DIM> fine_mesh_generator;
203  fine_mesh_generator.SetDomain(extractor.GetOutput());
204  fine_mesh_generator.SetMaxElementArea(mElementSize);
205  fine_mesh_generator.SetHoles(holes);
206  fine_mesh_generator.SetRegionMarkers(regions);
207  if(mpDomain)
208  {
209  fine_mesh_generator.SetDomain(mpDomain);
210  }
211  fine_mesh_generator.Update();
212  boost::shared_ptr<DiscreteContinuumMesh<DIM, DIM> > p_fine_mesh = fine_mesh_generator.GetMesh();
213  mMesh = p_fine_mesh;
214  }
215  else
216  {
217  EXCEPTION("3d image meshing not yet implemented");
218  }
219 }
220 
221 // Explicit instantiation
222 template class ImageToMesh<2> ;
223 template class ImageToMesh<3> ;
boost::shared_ptr< DiscreteContinuumMesh< DIM, DIM > > mMesh
Definition: ImageToMesh.hpp:65
void SetTissueDomain(boost::shared_ptr< Part< DIM > > pTissueDomain)
Definition: ImageToMesh.cpp:90
boost::shared_ptr< Part< DIM > > mpDomain
Definition: ImageToMesh.hpp:70
void SetDomain(boost::shared_ptr< Part< SPACE_DIM > > pDomain)
static boost::shared_ptr< ImageToMesh< DIM > > Create()
Definition: ImageToMesh.cpp:64
void SetSmoothingLength(double value)
void SetInput(vtkSmartPointer< vtkImageData > pImage)
vtkSmartPointer< vtkImageData > mpImage
Definition: ImageToMesh.hpp:55
vtkSmartPointer< vtkPolyData > GetOutput()
boost::shared_ptr< DiscreteContinuumMesh< DIM, DIM > > GetMesh()
Definition: ImageToMesh.cpp:77
void SetRegionMarkers(std::vector< DimensionalChastePoint< SPACE_DIM > > regionMarkers)
void SetMaxElementArea(units::quantity< unit::volume > maxElementArea)
units::quantity< unit::volume > mElementSize
Definition: ImageToMesh.hpp:60
boost::shared_ptr< DiscreteContinuumMesh< ELEMENT_DIM, SPACE_DIM > > GetMesh()
void SetHoles(std::vector< DimensionalChastePoint< SPACE_DIM > > holes)
void SetDoSmoothing(bool doSmoothing)
void SetElementSize(units::quantity< unit::volume > elementSize)
Definition: ImageToMesh.cpp:96
void SetInputRaw(vtkImageData *pImage)
void SetInput(vtkSmartPointer< vtkPolyData > pInputSurface)
Definition: Part.hpp:60