Skip to content
Snippets Groups Projects
Forked from agnumpde / dune-fufem
720 commits behind the upstream repository.
  • Ansgar Burchardt's avatar
    c6360609
    VTKBasisGridFunction: do not check size in constructor · c6360609
    Ansgar Burchardt authored
    With a Dune::VTKSequenceWriter I would like to have code like
    
        Vector x;
        VTKSequenceWriter<...> writer(...);
        VTKBasisGridFunction<...> function(..., x, ...);
        writer.addVertexData(function);
    
        for (int i = 0; i < level; ++i) {
            ... // refine grid, solve PDE, change size of x
            writer.write(static_cast<double>(i));
        }
    
    As the size of "x" is only correct in the for-loop, it cannot be checked
    when constructing the "writer". However the "function" needs to be setup
    outside the loop as well. So the VTKBasisGridFunction constructor should
    not check the size of "x".
    c6360609
    History
    VTKBasisGridFunction: do not check size in constructor
    Ansgar Burchardt authored
    With a Dune::VTKSequenceWriter I would like to have code like
    
        Vector x;
        VTKSequenceWriter<...> writer(...);
        VTKBasisGridFunction<...> function(..., x, ...);
        writer.addVertexData(function);
    
        for (int i = 0; i < level; ++i) {
            ... // refine grid, solve PDE, change size of x
            writer.write(static_cast<double>(i));
        }
    
    As the size of "x" is only correct in the for-loop, it cannot be checked
    when constructing the "writer". However the "function" needs to be setup
    outside the loop as well. So the VTKBasisGridFunction constructor should
    not check the size of "x".
Code owners
Assign users and groups as approvers for specific file changes. Learn more.
vtkbasisgridfunction.hh 2.29 KiB
// -*- c-basic-offset:4 tab-width:4 -*-

#ifndef VTK_BASIS_GRID_FUNCTION_HH
#define VTK_BASIS_GRID_FUNCTION_HH

#include <dune/grid/io/file/vtk/function.hh>

#include <dune/fufem/functions/basisgridfunction.hh>

/** \brief A VTK basis grid function.
 *
 *  This function "evaluates" by evaluating the global basis and interpolating the function values.
 *  \tparam Basis   The global basis.
 *  \tparam CoefficientType The vector type for the coefficients.
*/
template<class Basis, class CoefficientType>
class VTKBasisGridFunction : public Dune::VTKFunction<typename Basis::GridView>
{
    typedef Dune::VTKFunction<typename Basis::GridView> Base;
    typedef typename  CoefficientType::value_type RangeType;
public:
    typedef typename Base::Entity Entity;
    typedef typename Base::ctype ctype;
    using Base::dim;

    /** \brief Construct from given global basis, coefficient vector and name.
     *
     *  \param basis    The global basis.
     *  \param v    A corresponding vector of coefficients.
     *  \param s    A name of the function.
     */
    VTKBasisGridFunction(const Basis &basis, const CoefficientType &v, const std::string &s) :
        basis_(basis),
        coeffs_(v),
        s_( s )
    {
        /* Nothing. */
    }

    /** \brief Get the number of components the function has. */
    virtual int ncomps () const
    {
        return CoefficientType::value_type::dimension;
    }

    /** \brief Locally evaluate a component of the function.
     *
     *  \param comp The component to evaluate.
     *  \param e    The element the local coordinates are taken from.
     *  \param xi   The local coordinates where to evaluate the function.
     */
    virtual double evaluate (int comp, const Entity &e,
                             const Dune::FieldVector<ctype,dim> &xi) const
    {
        RangeType y;

        // use LocalEvaluator to make function also work for nonconforming bases
        LocalEvaluator<Basis, CoefficientType, RangeType>::evaluateLocal(basis_, coeffs_, e, xi, y);

        return y[comp];
    }

    /** \brief Get the name of that function. */
    virtual std::string name () const
    {
        return s_;
    }

    /** \brief Destructor. */
    virtual ~VTKBasisGridFunction() {}

private:
    const Basis &basis_;
    const CoefficientType &coeffs_;
    std::string s_;
};
#endif