Chaste  Build::
VesselNetworkWriter.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 #ifdef CHASTE_VTK
37 #define _BACKWARD_BACKWARD_WARNING_H 1
38 #include <vtkDoubleArray.h>
39 #include <vtkCellData.h>
40 #include <vtkCellArray.h>
41 #include <vtkPointData.h>
42 #include <vtkLine.h>
43 #include <vtkXMLPolyDataWriter.h>
44 #include <vtkVersion.h>
45 #endif // CHASTE_VTK
46 #include "SmartPointers.hpp"
47 #include "Exception.hpp"
48 #include "PetscTools.hpp"
49 #include "VesselNetworkWriter.hpp"
50 
51 template <unsigned DIM>
53  mpVesselNetwork(),
54  mpVtkVesselNetwork(vtkSmartPointer<vtkPolyData>::New()),
55  mIsVtkNetworkUpToDate(false),
56  mFilename(),
57  mReferenceLength(1.e-6 * unit::metres)
58 {
59 
60 }
61 
62 template <unsigned DIM>
64 {
65 
66 }
67 
68 template <unsigned DIM>
69 boost::shared_ptr<VesselNetworkWriter<DIM> > VesselNetworkWriter<DIM>::Create()
70 {
71  MAKE_PTR(VesselNetworkWriter<DIM>, pSelf);
72  return pSelf;
73 }
74 
75 template <unsigned DIM>
77 {
78  mpVesselNetwork = pNetwork;
79  mIsVtkNetworkUpToDate = false;
80 }
81 
82 template <unsigned DIM>
83 void VesselNetworkWriter<DIM>::SetReferenceLengthScale(units::quantity<unit::length> rReferenceLength)
84 {
85  mReferenceLength = rReferenceLength;
86 }
87 
88 template <unsigned DIM>
89 vtkSmartPointer<vtkPolyData> VesselNetworkWriter<DIM>::GetOutput()
90 {
91  if(!mpVesselNetwork)
92  {
93  EXCEPTION("A vessel network is required for the vtk writer.");
94  }
95 
96  if(mReferenceLength == 0.0 * unit::metres)
97  {
98  EXCEPTION("A non zero reference length scale is required for the vtk writer.");
99  }
100 
101  mpVtkVesselNetwork = vtkSmartPointer<vtkPolyData>::New();
102 
103  if(mpVesselNetwork->GetNumberOfVessels()>0)
104  {
105  // Set up the vessel data arrays.
106  std::vector<boost::shared_ptr<Vessel<DIM> > > vessels = mpVesselNetwork->GetVessels();
107  std::vector<vtkSmartPointer<vtkDoubleArray> > pVesselInfoVector;
108  std::map<std::string, double>::iterator vessel_map_iterator;
109  std::map<std::string, double> vessel_data_map = vessels[0]->GetOutputData();
110  for(vessel_map_iterator = vessel_data_map.begin(); vessel_map_iterator != vessel_data_map.end(); vessel_map_iterator++)
111  {
112  vtkSmartPointer<vtkDoubleArray> pVesselInfo = vtkSmartPointer<vtkDoubleArray>::New();
113  pVesselInfo->SetNumberOfComponents(1);
114  pVesselInfo->SetNumberOfTuples(vessels.size());
115  pVesselInfo->SetName((*vessel_map_iterator).first.c_str());
116  pVesselInfoVector.push_back(pVesselInfo);
117  }
118 
119  // Set up node data arrays
120  std::vector<boost::shared_ptr<VesselNode<DIM> > > nodes = mpVesselNetwork->GetNodes();
121  std::vector<vtkSmartPointer<vtkDoubleArray> > pNodeInfoVector;
122  std::map<std::string, double>::iterator vtk_node_map_iterator;
123  std::map<std::string, double> vtk_node_map = vessels[0]->GetStartNode()->GetOutputData();
124  for(vtk_node_map_iterator = vtk_node_map.begin(); vtk_node_map_iterator != vtk_node_map.end(); vtk_node_map_iterator++)
125  {
126  vtkSmartPointer<vtkDoubleArray> pNodeInfo = vtkSmartPointer<vtkDoubleArray>::New();
127  pNodeInfo->SetNumberOfComponents(1);
128  pNodeInfo->SetNumberOfTuples(nodes.size());
129  pNodeInfo->SetName((*vtk_node_map_iterator).first.c_str());
130  pNodeInfoVector.push_back(pNodeInfo);
131  }
132 
133  // Create the geometric data
134  vtkSmartPointer<vtkPoints> pPoints= vtkSmartPointer<vtkPoints>::New();
135  vtkSmartPointer<vtkCellArray> pLines = vtkSmartPointer<vtkCellArray>::New();
136  unsigned vessel_index=0;
137  for(unsigned idx=0; idx<nodes.size(); idx++)
138  {
139  nodes[idx]->SetId(idx);
140  if(DIM == 2)
141  {
142  pPoints->InsertNextPoint(nodes[idx]->rGetLocation().GetLocation(mReferenceLength)[0],
143  nodes[idx]->rGetLocation().GetLocation(mReferenceLength)[1], 0.0);
144  }
145  else
146  {
147  pPoints->InsertNextPoint(nodes[idx]->rGetLocation().GetLocation(mReferenceLength)[0],
148  nodes[idx]->rGetLocation().GetLocation(mReferenceLength)[1],
149  nodes[idx]->rGetLocation().GetLocation(mReferenceLength)[2]);
150  }
151 
152  std::map<std::string, double> vtk_node_data = nodes[idx]->GetOutputData();
153  for(unsigned jdx=0; jdx < pNodeInfoVector.size(); jdx++)
154  {
155  std::string key = pNodeInfoVector[jdx]->GetName();
156  // check if the key exists in the data map, if yes populate the vtk node data array
157  if(vtk_node_data.count(key) == 1)
158  {
159  pNodeInfoVector[jdx]->SetValue(idx, vtk_node_data[key]);
160  }
161  }
162  }
163 
164  // Create the vessels
165  typename std::vector<boost::shared_ptr<Vessel<DIM> > >::iterator it;
166  for(it = vessels.begin(); it < vessels.end(); it++)
167  {
168  vtkSmartPointer<vtkLine> pLine = vtkSmartPointer<vtkLine>::New();
169  std::vector<boost::shared_ptr<VesselSegment<DIM> > > segments = (*it)->GetSegments();
170  for(unsigned i = 0; i < segments.size(); i++)
171  {
172  pLine->GetPointIds()->InsertId(i, segments[i]->GetNode(0)->GetId());
173 
174  // Do an extra insert for the last node in the segment
175  if (i == segments.size() - 1)
176  {
177  pLine->GetPointIds()->InsertId(i + 1, segments[i]->GetNode(1)->GetId());
178  }
179  }
180  pLines->InsertNextCell(pLine);
181 
182  // Add the vessel data
183  std::map<std::string, double> vtk_vessel_data = (*it)->GetOutputData();
184  for(unsigned idx=0; idx < pVesselInfoVector.size(); idx++)
185  {
186  std::string key = pVesselInfoVector[idx]->GetName();
187  // If the key is in the vtk data use it
188  if(vtk_vessel_data.count(key) == 1)
189  {
190  pVesselInfoVector[idx]->SetValue(vessel_index, vtk_vessel_data[key]);
191  }
192  }
193  vessel_index++;
194  }
195  mpVtkVesselNetwork->SetPoints(pPoints);
196  mpVtkVesselNetwork->SetLines(pLines);
197 
198  for (unsigned i = 0; i < pVesselInfoVector.size(); i++)
199  {
200  mpVtkVesselNetwork->GetCellData()->AddArray(pVesselInfoVector[i]);
201  }
202 
203  for (unsigned i = 0; i < pNodeInfoVector.size(); i++)
204  {
205  mpVtkVesselNetwork->GetPointData()->AddArray(pNodeInfoVector[i]);
206  }
207  }
208  mIsVtkNetworkUpToDate = true;
209  return mpVtkVesselNetwork;
210 }
211 
212 template <unsigned DIM>
213 void VesselNetworkWriter<DIM>::SetFileName(const std::string& rFileName)
214 {
215  mFilename = rFileName;
216 }
217 
218 template <unsigned DIM>
220 {
221  if(mFilename.empty())
222  {
223  EXCEPTION("No file name set for VesselNetworkWriter");
224  }
225 
226  if(PetscTools::AmMaster())
227  {
228  vtkSmartPointer<vtkXMLPolyDataWriter> writer = vtkSmartPointer<vtkXMLPolyDataWriter>::New();
229  writer->SetFileName(mFilename.c_str());
230 
231  #if VTK_MAJOR_VERSION <= 5
232  writer->SetInput(this->GetOutput());
233  #else
234  writer->SetInputData(this->GetOutput());
235  #endif
236  writer->Write();
237  }
238 }
239 
240 // Explicit instantiation
241 template class VesselNetworkWriter<2>;
242 template class VesselNetworkWriter<3>;
boost::shared_ptr< VesselNetwork< DIM > > mpVesselNetwork
static boost::shared_ptr< VesselNetworkWriter< DIM > > Create()
void SetVesselNetwork(boost::shared_ptr< VesselNetwork< DIM > > pNetwork)
void SetReferenceLengthScale(units::quantity< unit::length > rReferenceLength)
vtkSmartPointer< vtkPolyData > mpVtkVesselNetwork
void SetFileName(const std::string &rFileName)
boost::shared_ptr< VesselNode< DIM > > GetNode(unsigned index)
Definition: Vessel.cpp:574
units::quantity< unit::length > mReferenceLength
vtkSmartPointer< vtkPolyData > GetOutput()