Skip to content
Snippets Groups Projects
Commit f59b1c41 authored by Peter Bastian's avatar Peter Bastian
Browse files

added VTK output to functions

[[Imported from SVN: r4391]]
parent 0af86f5b
Branches
Tags
No related merge requests found
......@@ -10,6 +10,8 @@
#include "common/exceptions.hh"
#include "grid/common/grid.hh"
#include "grid/utility/hierarchicsearch.hh"
#include "dune/io/file/vtk/vtkwriter.hh"
/**
* @file dune/disc/functions/functions.hh
......@@ -518,6 +520,53 @@ namespace Dune
};
// Wrapper for VTK output of grid functions
template<class GridImp, class IS, class RT, int m>
class VTKGridFunctionWrapper : public VTKWriter<GridImp,IS>::VTKFunction
{
enum {n=GridImp::dimension};
typedef typename GridImp::ctype DT;
typedef typename GridImp::Traits::template Codim<0>::Entity Entity;
public:
//! return number of components
virtual int ncomps () const
{
return m;
}
//! evaluate single component comp in the entity e at local coordinates xi
/*! Evaluate the function in an entity at local coordinates.
@param[in] comp number of component to be evaluated
@param[in] e reference to grid entity of codimension 0
@param[in] xi point in local coordinates of the reference element of e
\return value of the component
*/
virtual double evaluate (int comp, const Entity& e, const Dune::FieldVector<DT,n>& xi) const
{
return func.evallocal(comp,e,xi);
}
// get name
virtual std::string name () const
{
return myname;
}
// constructor remembers reference to a grid function
VTKGridFunctionWrapper (const GridFunction<GridImp,RT,m>& f, std::string s) : func(f), myname(s)
{}
virtual ~VTKGridFunctionWrapper() {}
private:
const GridFunction<GridImp,RT,m>& func;
std::string myname;
};
/** @} */
......
......@@ -22,6 +22,7 @@
#include "dune/istl/operators.hh"
#include "dune/istl/bcrsmatrix.hh"
#include "dune/disc/shapefunctions/lagrangeshapefunctions.hh"
#include "dune/io/file/vtk/vtkwriter.hh"
// same directory includes
#include "functions.hh"
......@@ -97,7 +98,7 @@ namespace Dune
std::cerr << "not enough memory in P0Function" << std::endl;
throw; // rethrow exception
}
std::cout << "making function with " << mapper_.size() << " components" << std::endl;
// std::cout << "making function with " << mapper_.size() << " components" << std::endl;
}
//! create function with initialization from a vector
......@@ -113,7 +114,7 @@ namespace Dune
std::cerr << "not enough memory in P0Function" << std::endl;
throw; // rethrow exception
}
std::cout << "making function with " << mapper_.size() << " components" << std::endl;
// std::cout << "making function with " << mapper_.size() << " components" << std::endl;
// check size of initialization vector
if (v.size()!=(unsigned int)mapper_.size())
......@@ -251,6 +252,13 @@ namespace Dune
return mapper_;
}
//! VTK output
void vtkout (VTKWriter<G,IS>& vtkwriter, std::string s)
{
typename VTKWriter<G,IS>::VTKFunction *p = new VTKGridFunctionWrapper<G,IS,RT,m>(*this,s);
vtkwriter.addCellData(p);
}
private:
// a reference to the grid
const G& grid_;
......
......@@ -26,6 +26,7 @@
#include "dune/istl/bcrsmatrix.hh"
#include "dune/istl/owneroverlapcopy.hh"
#include "dune/disc/shapefunctions/lagrangeshapefunctions.hh"
#include "dune/io/file/vtk/vtkwriter.hh"
// same directory includes
#include "functions.hh"
......@@ -929,6 +930,13 @@ namespace Dune
return mapper_;
}
//! VTK output
void vtkout (VTKWriter<G,IS>& vtkwriter, std::string s)
{
typename VTKWriter<G,IS>::VTKFunction *p = new VTKGridFunctionWrapper<G,IS,RT,m>(*this,s);
vtkwriter.addVertexData(p);
}
private:
// a reference to the grid
const G& grid_;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment