Skip to content
Snippets Groups Projects
Commit 622a4862 authored by Dr. Felix Tobias Schindler's avatar Dr. Felix Tobias Schindler Committed by Tobias Leibner
Browse files

[functions] add visualize()

parent 36950d66
No related branches found
No related tags found
1 merge request!44ccache and clang sanitizer flags for CI
......@@ -202,22 +202,24 @@ private:
using DomainType = typename LocalFunctionType::DomainType;
public:
VisualizationAdapter(const GridFunctionType& localizable_function,
VisualizationAdapter(const GridFunctionType& grid_function,
const VisualizerInterface<range_dim, range_dim_cols, RangeField>& visualizer,
const std::string nm = "",
const XT::Common::Parameter& param = {})
: local_function_(localizable_function.local_function())
: function_(grid_function.copy_as_grid_function())
, local_function_(function_->local_function())
, visualizer_(visualizer)
, name_(nm.empty() ? localizable_function.name() : nm)
, name_(nm.empty() ? function_->name() : nm)
, param_(param)
{}
VisualizationAdapter(const GridFunctionType& localizable_function,
VisualizationAdapter(const GridFunctionType& grid_function,
const std::string nm = "",
const XT::Common::Parameter& param = {})
: local_function_(localizable_function.local_function())
: function_(grid_function.copy_as_grid_function())
, local_function_(function_->local_function())
, visualizer_(new DefaultVisualizer<range_dim, range_dim_cols, RangeField>())
, name_(nm.empty() ? localizable_function.name() : nm)
, name_(nm.empty() ? function_->name() : nm)
, param_(param)
{}
......@@ -239,12 +241,14 @@ public:
}
private:
const std::unique_ptr<GridFunctionType> function_;
mutable std::unique_ptr<LocalFunctionType> local_function_;
const Common::ConstStorageProvider<VisualizerInterface<range_dim, range_dim_cols, RangeField>> visualizer_;
const std::string name_;
const XT::Common::Parameter param_;
}; // class VisualizationAdapter
template <class GridViewType, size_t range_dim, size_t range_dim_cols, class RangeField>
class GradientVisualizationAdapter : public VTKFunction<GridViewType>
{
......@@ -261,25 +265,27 @@ private:
using DomainType = typename LocalFunctionType::DomainType;
public:
GradientVisualizationAdapter(const GridFunctionType& localizable_function,
GradientVisualizationAdapter(const GridFunctionType& grid_function,
const VisualizerInterface<d, 1, RangeField>& visualizer,
const std::string nm = "",
const XT::Common::Parameter& param = {})
: local_function_(localizable_function.local_function())
: function_(grid_function.copy_as_grid_function())
, local_function_(function_->local_function())
, visualizer_(visualizer)
, name_(nm.empty() ? "grad_ " + localizable_function.name() : nm)
, name_(nm.empty() ? "grad_ " + function_->name() : nm)
, param_(param)
{
if (range_dim > 1)
DUNE_THROW(Dune::NotImplemented, "Only implemented for scalar functions by now!");
}
GradientVisualizationAdapter(const GridFunctionType& localizable_function,
GradientVisualizationAdapter(const GridFunctionType& grid_function,
const std::string nm = "",
const XT::Common::Parameter& param = {})
: local_function_(localizable_function.local_function())
: function_(grid_function.copy_as_grid_function())
, local_function_(function_->local_function())
, visualizer_(new DefaultVisualizer<d, 1, RangeField>())
, name_(nm.empty() ? "grad_" + localizable_function.name() : nm)
, name_(nm.empty() ? "grad_" + function_->name() : nm)
, param_(param)
{
if (range_dim > 1)
......@@ -304,6 +310,7 @@ public:
}
private:
const std::unique_ptr<GridFunctionType> function_;
mutable std::unique_ptr<LocalFunctionType> local_function_;
const Common::ConstStorageProvider<VisualizerInterface<d, 1, RangeField>> visualizer_;
const std::string name_;
......
// This file is part of the dune-xt project:
// https://github.com/dune-community/dune-xt
// Copyright 2009-2020 dune-xt developers and contributors. All rights reserved.
// License: Dual licensed as BSD 2-Clause License (http://opensource.org/licenses/BSD-2-Clause)
// or GPL-2.0+ (http://opensource.org/licenses/gpl-license)
// with "runtime exception" (http://www.dune-project.org/license.html)
// Authors:
// Felix Schindler (2020)
#ifndef DUNE_XT_FUNCTIONS_VISUALIZATION_HH
#define DUNE_XT_FUNCTIONS_VISUALIZATION_HH
#include <dune/xt/grid/type_traits.hh>
#include <dune/xt/functions/base/visualization.hh>
#include <dune/xt/functions/interfaces/function.hh>
#include <dune/xt/functions/interfaces/grid-function.hh>
#include <dune/xt/functions/grid-function.hh>
namespace Dune {
namespace XT {
namespace Functions {
namespace internal {
template <class GridViewType>
std::unique_ptr<VTKWriter<GridViewType>>
create_vtkwriter(const GridViewType& grid_view, const bool subsampling = true, const int subsampling_level = 2)
{
static_assert(Grid::is_view<GridViewType>::value);
return subsampling ? std::make_unique<SubsamplingVTKWriter<GridViewType>>(grid_view, subsampling_level)
: std::make_unique<VTKWriter<GridViewType>>(grid_view, VTK::nonconforming);
}
template <class GridViewType, class E, size_t r, size_t rC, class R>
void add_to_vtkwriter(VTKWriter<GridViewType>& vtk_writer,
const GridFunctionInterface<E, r, rC, R>& grid_function,
const XT::Common::Parameter& param = {},
const VisualizerInterface<r, rC, R>& visualizer = default_visualizer<r, rC, R>())
{
static_assert(Grid::is_view<GridViewType>::value);
static_assert(std::is_same<E, XT::Grid::extract_entity_t<GridViewType>>::value);
vtk_writer.addVertexData(
std::make_shared<VisualizationAdapter<GridViewType, r, rC, R>>(grid_function, visualizer, "", param));
}
template <class GridViewType, size_t d, size_t r, size_t rC, class R>
void add_to_vtkwriter(VTKWriter<GridViewType>& vtk_writer,
const FunctionInterface<d, r, rC, R>& function,
const XT::Common::Parameter& param = {},
const VisualizerInterface<r, rC, R>& visualizer = default_visualizer<r, rC, R>())
{
static_assert(Grid::is_view<GridViewType>::value);
static_assert(GridViewType::dimension == d);
add_to_vtkwriter(vtk_writer, make_grid_function<GridViewType>(function), param, visualizer);
}
template <class GridViewType, class E, size_t r, size_t rC, class R>
void add_gradient_to_vtkwriter(VTKWriter<GridViewType>& vtk_writer,
const GridFunctionInterface<E, r, rC, R>& grid_function,
const XT::Common::Parameter& param = {},
const VisualizerInterface<GridViewType::dimension, 1, R>& visualizer =
default_visualizer<GridViewType::dimension, 1, R>())
{
static_assert(Grid::is_view<GridViewType>::value);
static_assert(std::is_same<E, XT::Grid::extract_entity_t<GridViewType>>::value);
const auto adapter =
std::make_shared<GradientVisualizationAdapter<GridViewType, r, rC, R>>(grid_function, visualizer, "", param);
vtk_writer.addCellData(adapter);
}
template <class GridViewType, size_t d, size_t r, size_t rC, class R>
void add_gradient_to_vtkwriter(VTKWriter<GridViewType>& vtk_writer,
const FunctionInterface<d, r, rC, R>& function,
const XT::Common::Parameter& param = {},
const VisualizerInterface<GridViewType::dimension, 1, R>& visualizer =
default_visualizer<GridViewType::dimension, 1, R>())
{
static_assert(Grid::is_view<GridViewType>::value);
static_assert(GridViewType::dimension == d);
add_gradient_to_vtkwriter(vtk_writer, make_grid_function<GridViewType>(function), param, visualizer);
}
template <class GridViewType>
auto write_visualization(VTKWriter<GridViewType>& vtk_writer,
const std::string path,
const VTK::OutputType vtk_output_type = VTK::appendedraw)
{
static_assert(Grid::is_view<GridViewType>::value);
if (path.empty())
DUNE_THROW(Exceptions::wrong_input_given, "path must not be empty!");
const auto directory = Common::directory_only(path);
Common::test_create_directory(directory);
if (MPIHelper::getCollectiveCommunication().size() == 1)
vtk_writer.write(path, vtk_output_type);
else
vtk_writer.pwrite(Common::filename_only(path), directory, "", vtk_output_type);
}
} // namespace internal
/**
* \note We use the SubsamplingVTKWriter (which is better for higher orders) by default: the grid you see in the
* visualization may thus be a refinement of the actual grid!
*/
template <class GridViewType, size_t d, size_t r, size_t rC, class R>
void visualize(const FunctionInterface<d, r, rC, R>& function,
const GridViewType& grid_view,
const std::string path,
const bool subsampling = true,
const VTK::OutputType vtk_output_type = VTK::appendedraw,
const XT::Common::Parameter& param = {},
const VisualizerInterface<r, rC, R>& visualizer = default_visualizer<r, rC, R>())
{
static_assert(Grid::is_view<GridViewType>::value);
static_assert(GridViewType::dimension == d);
auto vtk_writer = internal::create_vtkwriter(grid_view, subsampling);
internal::add_to_vtkwriter(*vtk_writer, function, param, visualizer);
internal::write_visualization(*vtk_writer, path, vtk_output_type);
} // ... visualize(...)
/**
* \note We use the SubsamplingVTKWriter (which is better for higher orders) by default: the grid you see in the
* visualization may thus be a refinement of the actual grid!
*/
template <class GridViewType, class E, size_t r, size_t rC, class R>
void visualize(const GridFunctionInterface<E, r, rC, R>& grid_function,
const GridViewType& grid_view,
const std::string path,
const bool subsampling = true,
const VTK::OutputType vtk_output_type = VTK::appendedraw,
const XT::Common::Parameter& param = {},
const VisualizerInterface<r, rC, R>& visualizer = default_visualizer<r, rC, R>())
{
static_assert(Grid::is_view<GridViewType>::value);
static_assert(std::is_same<E, XT::Grid::extract_entity_t<GridViewType>>::value);
auto vtk_writer = internal::create_vtkwriter(grid_view, subsampling);
internal::add_to_vtkwriter(*vtk_writer, grid_function, param, visualizer);
internal::write_visualization(*vtk_writer, path, vtk_output_type);
} // ... visualize(...)
/**
* \note We use the SubsamplingVTKWriter (which is better for higher orders) by default: the grid you see in the
* visualization may thus be a refinement of the actual grid!
* \note Not yet implemented for vector-valued functions.
*/
template <class GridViewType, size_t d, size_t r, size_t rC, class R>
void visualize_gradient(const FunctionInterface<d, r, rC, R>& function,
const GridViewType& grid_view,
const std::string path,
const bool subsampling = true,
const VTK::OutputType vtk_output_type = VTK::appendedraw,
const XT::Common::Parameter& param = {},
const VisualizerInterface<d, 1, R>& visualizer = default_visualizer<d, 1, R>())
{
static_assert(Grid::is_view<GridViewType>::value);
auto vtk_writer = internal::create_vtkwriter(grid_view, subsampling);
add_gradient_to_vtkwriter(*vtk_writer, function, param, visualizer);
internal::write_visualization(*vtk_writer, path, vtk_output_type);
} // ... visualize_gradient(...)
/**
* \note We use the SubsamplingVTKWriter (which is better for higher orders) by default: the grid you see in the
* visualization may thus be a refinement of the actual grid!
* \note Not yet implemented for vector-valued functions.
*/
template <class GridViewType, class E, size_t r, size_t rC, class R>
void visualize_gradient(const GridFunctionInterface<E, r, rC, R>& grid_function,
const GridViewType& grid_view,
const std::string path,
const bool subsampling = true,
const VTK::OutputType vtk_output_type = VTK::appendedraw,
const XT::Common::Parameter& param = {},
const VisualizerInterface<GridViewType::dimension, 1, R>& visualizer =
default_visualizer<GridViewType::dimension, 1, R>())
{
static_assert(Grid::is_view<GridViewType>::value);
auto vtk_writer = internal::create_vtkwriter(grid_view, subsampling);
add_gradient_to_vtkwriter(*vtk_writer, grid_function, param, visualizer);
internal::write_visualization(*vtk_writer, path, vtk_output_type);
} // ... visualize_gradient(...)
} // namespace Functions
} // namespace XT
} // namespace Dune
#endif // DUNE_XT_FUNCTIONS_VISUALIZATION_HH
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment