Mercurial > hg > kwantix
changeset 42:3f8311cbf602
Setup initial class for vtkplot; cleanup interpolator to use C++1x's auto and untabify
author | Jordi Gutiérrez Hermoso <jordigh@gmail.com> |
---|---|
date | Sat, 13 Mar 2010 15:55:11 -0600 |
parents | 834f8e778a65 |
children | 8aea477c7cf8 |
files | src/include/interpolator.hpp src/include/vtkplot.hpp src/interpolator.cpp src/vtkplot.cpp |
diffstat | 4 files changed, 88 insertions(+), 64 deletions(-) [+] |
line wrap: on
line diff
--- a/src/include/interpolator.hpp +++ b/src/include/interpolator.hpp @@ -25,6 +25,8 @@ using std::map; using boost::shared_ptr; +class vtkplot; + //FIXME: Break this up into an ansatz class and an interpolator //class. /*! \brief An interpolator \f$ \displaystyle u(x) = \sum_{\xi \in @@ -190,6 +192,9 @@ * values that the interpolator takes on that point. */ void into_os(std::ostream& os) const; + + ///Needed in order to plot the data + friend class vtkplot; private: ///Once the matrix is defined, this function inverts it.
--- a/src/include/vtkplot.hpp +++ b/src/include/vtkplot.hpp @@ -10,7 +10,6 @@ namespace kwantix{ ///A class for creating VTK plots -template<typename RBF> class vtkplot{ public: /*! \brief Constructor performs a Delaunay triangulation on data @@ -19,11 +18,13 @@ * Delaunay triangulation will be done, and the * third column is the value at this point. */ + template<typename RBF> vtkplot(const interpolator<RBF>& u, - bool offscreen = true); + bool offscreen_in = false); ///Defines the view direction when plotting. void set_view_direction(const vector& dir); ///Without changing the triangulation, update the values. + template<typename RBF> void update_values(const interpolator<RBF>& u); ///Enable offscreen plotting. void plot_offscreen(bool offscreen = true); @@ -33,6 +34,19 @@ void plot() const; private: + + ///Setup the initial pipeline + template<typename RBF> + void SetupPipeline(const interpolator<RBF>& u); + + ///Setup the pipeline to plot offscreen or not + void InitOffscreen(); + + template<typename RBF> + static void AdjustPoints(void* arguments); + + ///The associated interpolator's hash + size_t hash; ///Base filename string basename; ///Direction from which to plot
--- a/src/interpolator.cpp +++ b/src/interpolator.cpp @@ -1,9 +1,5 @@ #include <vector> #include <set> - -//debug -#include <iostream> - #include <utility> #include "include/interpolator.hpp" @@ -104,17 +100,20 @@ shared_ptr<dirichlet_op> D(new dirichlet_op); set<point> intr; set<point> bdry; //empty - map<point, point> nrml; //empty + map<point, vector> nrml; //empty map<point, double> g; //empty bool dim_set = false; size_t dimension = 0; - for(map<point,double>::const_iterator i = Xi.begin(); i != Xi.end(); i++){ - if(!dim_set){ + for(auto i = Xi.begin(); i != Xi.end(); i++) + { + if(!dim_set) + { dimension = (i -> first).size(); dim_set = true; } - else if(dimension != (i -> first).size()){ + else if(dimension != (i -> first).size()) + { badArgument exc; exc.reason = "Inconformant dimensions in interpolation data."; exc.line = __LINE__;
--- a/src/vtkplot.cpp +++ b/src/vtkplot.cpp @@ -44,21 +44,11 @@ #include <sstream> #include <iomanip> -// *************** move this stuff into a class ************************ - -unsigned int GridSize = 50; - -using namespace std; - -double F(double x, double y, double t) -{ - return sin(y)*cos(x)*sin(t); -} +namespace kwantix{ class CommandSubclass : public vtkCommand { public: - vtkTypeRevisionMacro(CommandSubclass, vtkCommand); static CommandSubclass *New() { return new CommandSubclass; @@ -67,8 +57,6 @@ void Execute(vtkObject *caller, unsigned long eventId, void *callData) { - cout << "timer callback" << endl; - vtkRenderWindowInteractor *iren = static_cast<vtkRenderWindowInteractor*>(caller); @@ -82,14 +70,15 @@ }; -vtkCxxRevisionMacro(CommandSubclass, "$Revision: 1.1 $"); - -void AdjustPoints(void* arguments) +template<typename RBF> +void vtkplot::AdjustPoints(void* arguments) { - cout << "AdjustPoints" << endl; vtkProgrammableFilter* programmableFilter = static_cast<vtkProgrammableFilter*>(arguments); - + arguments = programmableFilter + 1; + + auto u = static_cast< const interpolator<RBF>& > (arguments); + vtkPoints* inPts = programmableFilter->GetPolyDataInput()->GetPoints(); vtkIdType numPts = inPts->GetNumberOfPoints(); vtkSmartPointer<vtkPoints> newPts = @@ -101,11 +90,13 @@ static double t = 0; + //auto& vals = u.at(); + for(vtkIdType i = 0; i < numPts; i++) { double p[3]; inPts -> GetPoint(i,p); - p[2] = F(p[0],p[1],t); + //p[2] = F(p[0],p[1],t); newPts -> SetPoint(i,p); newScalars -> InsertTuple(i, &p[2]); } @@ -121,46 +112,43 @@ ->GetPointData() -> SetScalars(newScalars); } -typedef double (*Func3d)(double, double, double); - -void InitOffscreen(bool offscreen); - -vtkSmartPointer<vtkDelaunay2D> TriangulateTerrain(Func3d func); - -void PlotPoints(vtkSmartPointer<vtkDelaunay2D> delaunay, - bool offscreen); - -void InitOffscreen(bool offscreen) +template<typename RBF> +vtkplot::vtkplot(const interpolator<RBF>& u, + bool offscreen_in) + :hash(u.rbfs_hash), offscreen(offscreen) { - // Graphics Factory - vtkSmartPointer<vtkGraphicsFactory> graphics_factory - = vtkSmartPointer<vtkGraphicsFactory>::New(); - graphics_factory->SetOffScreenOnlyMode( offscreen ); - graphics_factory->SetUseMesaClasses( offscreen ); - - // Imaging Factory - vtkSmartPointer<vtkImagingFactory> imaging_factory - = vtkSmartPointer<vtkImagingFactory>::New(); - imaging_factory->SetUseMesaClasses( offscreen ); + if(u.thebvp -> get_domain() -> get_dimension() != 2) + { + badArgument exc; + exc.reason = "This class only works on interpolators whose domain lies in R^2"; + exc.line = __LINE__; + exc.file = __FILE__; + throw exc; + } + SetupPipeline(u, offscreen); } -vtkSmartPointer<vtkDelaunay2D> TriangulateTerrain(Func3d func) +template<typename RBF> +void vtkplot::SetupPipeline(const interpolator<RBF>& u) { vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New(); vtkSmartPointer<vtkDoubleArray> scalars = vtkSmartPointer<vtkDoubleArray>::New(); - - for(double x = -3; x < 3; x += 6.0/GridSize) - { - for(double y = -3; y < 3; y += 6.0/GridSize) - { - double r = sqrt(9-y*y)*sqrt(9-x*x)/3; - double z = func(x,r,0); - vtkIdType idx = points->InsertNextPoint(x,r, z); - scalars -> InsertTuple(idx, &z); - } - } + for(auto i = u.thebvp -> get_domain() -> get_interior().begin(); + i != u.thebvp -> get_domain() -> get_interior().end(); i++) + { + double z = 0; + vtkIdType idx = points -> InsertNextPoint( *i(1), *i(2), 0); + scalars -> InsertTuple(idx, &z); + } + for(auto i = u.thebvp -> get_domain() -> get_boundary().begin(); + i != u.thebvp -> get_domain() -> get_boundary().end(); i++) + { + double z = 0; + vtkIdType idx = points -> InsertNextPoint( *i(1), *i(2), 0); + scalars -> InsertTuple(idx, &z); + } //add the grid points to a polydata object vtkSmartPointer<vtkPolyData> polydata = vtkSmartPointer<vtkPolyData>::New(); @@ -173,13 +161,12 @@ delaunay->SetInput(polydata); delaunay->Update(); - bool offscreen = false; - InitOffscreen(offscreen); + InitOffscreen(); vtkSmartPointer<vtkProgrammableFilter> programmableFilter = vtkSmartPointer<vtkProgrammableFilter>::New(); programmableFilter->SetInput(delaunay -> GetOutput()); - programmableFilter->SetExecuteMethod(AdjustPoints, programmableFilter); + programmableFilter->SetExecuteMethod(AdjustPoints<RBF>, programmableFilter, u); //Normals for Gourad shading @@ -252,6 +239,7 @@ { for(size_t i = 0; i < 1000; i++) { + using namespace std; string filename = "foo"; stringstream ss; string i_str; ss << setw(4) << setfill('0') << i; @@ -287,3 +275,21 @@ } return delaunay; } + +void vtkplot::InitOffscreen() +{ + // Graphics Factory + vtkSmartPointer<vtkGraphicsFactory> graphics_factory + = vtkSmartPointer<vtkGraphicsFactory>::New(); + graphics_factory->SetOffScreenOnlyMode( offscreen ); + graphics_factory->SetUseMesaClasses( offscreen ); + + // Imaging Factory + vtkSmartPointer<vtkImagingFactory> imaging_factory + = vtkSmartPointer<vtkImagingFactory>::New(); + imaging_factory->SetUseMesaClasses( offscreen ); +} + + +} //namespace kwantix +