Chaste  Build::
GeometryWriter.cpp
1 /*
2 
3 Copyright (c) 2005-2015, 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 #define _BACKWARD_BACKWARD_WARNING_H 1
37 #include <vtkDoubleArray.h>
38 #include <vtkCellData.h>
39 #include <vtkCellArray.h>
40 #include <vtkPointData.h>
41 #include <vtkLine.h>
42 #include <vtkXMLPolyDataWriter.h>
43 #include <vtkSTLWriter.h>
44 #include <vtkVersion.h>
45 #include <vtkGeometryFilter.h>
46 #include <vtkTriangleFilter.h>
47 #include <vtkCleanPolyData.h>
48 #include "SmartPointers.hpp"
49 #include "Exception.hpp"
50 #include "PetscTools.hpp"
51 #include "GeometryWriter.hpp"
52 
54  mpInputGeometry(),
55  mFilename(),
56  mFormat(GeometryFormat::VTP)
57 {
58 
59 }
60 
62 {
63 
64 }
65 
66 boost::shared_ptr<GeometryWriter > GeometryWriter::Create()
67 {
68  MAKE_PTR(GeometryWriter, pSelf);
69  return pSelf;
70 }
71 
72 void GeometryWriter::SetFileName(const std::string& rFileName)
73 {
74  mFilename = rFileName;
75 }
76 
77 void GeometryWriter::SetInput(vtkSmartPointer<vtkPolyData> pSurface)
78 {
79  mpInputGeometry = pSurface;
80 }
81 
83 {
84  mFormat = format;
85 }
86 
88 {
89  if(!mpInputGeometry)
90  {
91  EXCEPTION("An input geometry is not set.");
92  }
93 
94  if(mFilename.empty())
95  {
96  EXCEPTION("No file name set for the GeometryWriter.");
97  }
98 
99  if(mFormat == GeometryFormat::STL)
100  {
101  if(PetscTools::AmMaster())
102  {
103  vtkSmartPointer<vtkTriangleFilter> p_tri_filter = vtkSmartPointer<vtkTriangleFilter>::New();
104  #if VTK_MAJOR_VERSION <= 5
105  p_tri_filter->SetInput(mpInputGeometry);
106  #else
107  p_tri_filter->SetInputData(mpInputGeometry);
108  #endif
109 
110  vtkSmartPointer<vtkCleanPolyData> p_clean_filter = vtkSmartPointer<vtkCleanPolyData>::New();
111  p_clean_filter->SetInputConnection(p_tri_filter->GetOutputPort());
112  p_clean_filter->Update();
113 
114  vtkSmartPointer<vtkSTLWriter> writer = vtkSmartPointer<vtkSTLWriter>::New();
115  writer->SetFileName(mFilename.c_str());
116  #if VTK_MAJOR_VERSION <= 5
117  writer->SetInput(p_clean_filter->GetOutput());
118  #else
119  writer->SetInputData(p_clean_filter->GetOutput());
120  #endif
121  writer->SetFileTypeToASCII();
122  writer->Write();
123  }
124  }
125  else
126  {
127  if(PetscTools::AmMaster())
128  {
129  vtkSmartPointer<vtkXMLPolyDataWriter> writer = vtkSmartPointer<vtkXMLPolyDataWriter>::New();
130  writer->SetFileName(mFilename.c_str());
131  #if VTK_MAJOR_VERSION <= 5
132  writer->SetInput(mpInputGeometry);
133  #else
134  writer->SetInputData(mpInputGeometry);
135  #endif
136  writer->Write();
137  }
138  }
139 }
GeometryFormat::Value mFormat
static boost::shared_ptr< GeometryWriter > Create()
void SetOutputFormat(GeometryFormat::Value format)
vtkSmartPointer< vtkPolyData > mpInputGeometry
void SetInput(vtkSmartPointer< vtkPolyData > pSurface)
std::string mFilename
void SetFileName(const std::string &rFileName)