#ifndef DUNE_DORIE_GRID_FUNCTION_WRITER_HH #define DUNE_DORIE_GRID_FUNCTION_WRITER_HH #include #include #include #include #include #include #include #include namespace Dune{ namespace Dorie{ /** * @brief Converts PDELab grid function to match a dune-function. * @note This adapter differs from the one provided in PDELab in that this * one has the current interface of dune-grid to write vtk files. * The main difference using this adapter is that it allows to write * more than 3 components. * * @tparam GF PDELab grid function type. */ template class VTKGridFunctionAdapter { using GridFunction = typename std::decay::type; using GridView = typename GF::Traits::GridViewType; using Entity = typename GridView::template Codim<0>::Entity; using Range = typename GF::Traits::RangeType; using Domain = typename GF::Traits::DomainType; public: /** * @brief Constructs the object * * @param[in] gf pointer to a constant grid function */ VTKGridFunctionAdapter(std::shared_ptr gf) : _gf(gf) {} /** * @brief Binds an entity to the grid function * * @param[in] e Entity */ void bind(const Entity& e) { _entity = &e; } /** * @brief Unbinds the previously binded entity */ void unbind() { _entity = nullptr; } /** * @brief Evaluates the grid function at a given coordinate * * @param[in] x Local coordinate with respect to the last binded entity */ Range operator()(const Domain& x) const { Range y; _gf->evaluate(*_entity,x,y); return y; } private: const Entity* _entity; std::shared_ptr _gf; }; /*-------------------------------------------------------------------------*//** * @brief Class for grid fucntion vtk sequence writer. This class works * exacly the same way as a VTKSequence writer but it receives * directly smart pointers to grid functions * * @tparam GV GridView of the dune grid */ template class GridFunctionVTKSequenceWriter : public Dune::VTKSequenceWriter { public: // export constructors of vtk sequence writer using Dune::VTKSequenceWriter::VTKSequenceWriter; /*-----------------------------------------------------------------------*//** * @brief Adds a vertex data. * * @param[in] gf_ptr A pointer to a grid function of the type GF * @param[in] name Name of the variable * * @tparam GF Type of the grid function */ template void addVertexData(std::shared_ptr gf_ptr, std::string name) { static_assert(std::is_same::value, "GridFunctionVTKSequenceWriter only works with only one GridView type."); Dune::Dorie::VTKGridFunctionAdapter vtk_gf(gf_ptr); const int vtk_ncomps = GF::Traits::dimRange; const int dim = GF::Traits::dimDomain; auto vtk_type = (vtk_ncomps == dim) ? VTK::FieldInfo::Type::vector : VTK::FieldInfo::Type::scalar; auto vtk_info = VTK::FieldInfo(name,vtk_type,vtk_ncomps); this->vtkWriter()->addVertexData(vtk_gf,vtk_info); } /*-----------------------------------------------------------------------*//** * @brief Adds a cell data. * * @param[in] gf_ptr A pointer to a grid function of the type GF * @param[in] name Name of the variable * * @tparam GF Type of the grid function */ template void addCellData(std::shared_ptr gf_ptr, std::string name) { static_assert(std::is_same::value, "GridFunctionVTKSequenceWriter only works with only one GridView type."); Dune::Dorie::VTKGridFunctionAdapter vtk_gf(gf_ptr); const int vtk_ncomps = GF::Traits::dimRange; const int dim = GF::Traits::dimDomain; auto vtk_type = (vtk_ncomps == dim) ? VTK::FieldInfo::Type::vector : VTK::FieldInfo::Type::scalar; auto vtk_info = VTK::FieldInfo(name,vtk_type,vtk_ncomps); this->vtkWriter()->addCellData(vtk_gf,vtk_info); } /*-----------------------------------------------------------------------*//** * @brief Clear all the attached VTKFunctions. * @details This is necessary for each iteration if VTKFunctions are * managed internaly with pointers instead of references */ void clear() { this->vtkWriter()->clear(); } }; } // namespace Dorie } // namespace Dune #endif // DUNE_DORIE_GRID_FUNCTION_WRITER_HH