Skip to content

Commit

Permalink
Add VTK sample function and related test
Browse files Browse the repository at this point in the history
  • Loading branch information
pcafrica committed May 29, 2023
1 parent 048ea3e commit a202b0f
Show file tree
Hide file tree
Showing 5 changed files with 143 additions and 31 deletions.
46 changes: 15 additions & 31 deletions cmake/modules/FindDEAL_II_VTK.cmake
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
## ---------------------------------------------------------------------
##
## Copyright (C) 2022 by the deal.II authors
## Copyright (C) 2023 by the deal.II authors
##
## This file is part of the deal.II library.
##
Expand Down Expand Up @@ -39,43 +39,27 @@ if(VTK_FOUND)
set(VTK_MINOR_VERSION "${VTK_MINOR_VERSION}")

set(VTK_INCLUDE_DIR
${VTK_PREFIX_PATH}/include/vtk-${VTK_MAJOR_VERSION}.${VTK_MINOR_VERSION})
${VTK_PREFIX_PATH}/include/vtk-${VTK_MAJOR_VERSION}.${VTK_MINOR_VERSION})

#
# We'd like to have the full library names but the VTK package only
# exports a list with namespace-qualified names, e.g. VTK::CommonCore.
# So we check again for every lib and store the full path:
#
set(_libraries "")
# Try to find full paths from targets contained in VTK_LIBRARIES.
set(_libraries)
foreach(_library ${VTK_LIBRARIES})
if(NOT ${_library} MATCHES "Python"
AND NOT ${_library} MATCHES "MPI4Py")
get_target_property(_configurations ${_library} IMPORTED_CONFIGURATIONS)

foreach(_component ${VTK_AVAILABLE_COMPONENTS})
deal_ii_find_library(VTK_LIBRARY_${_component}
NAMES libvtk${_component}-${VTK_VERSION_MAJOR}.${VTK_VERSION_MINOR}.so
HINTS ${VTK_PREFIX_PATH}/lib
NO_DEFAULT_PATH
NO_CMAKE_ENVIRONMENT_PATH
NO_CMAKE_PATH
NO_SYSTEM_ENVIRONMENT_PATH
NO_CMAKE_SYSTEM_PATH
NO_CMAKE_FIND_ROOT_PATH
)

if (NOT "${VTK_LIBRARY_${_component}}" MATCHES "-NOTFOUND")
list(APPEND _libraries VTK_LIBRARY_${_component})
else()
# VTK_AVAILABLE_COMPONENTS contains also header-only modules.
# If the library has not been found, check if corresponding
# headers exist.
if (NOT EXISTS "${VTK_INCLUDE_DIR}/${_component}" AND
NOT EXISTS "${VTK_INCLUDE_DIR}/vtk_${_component}.h")
message(FATAL_ERROR "VTK: component ${_component} not found.")
if(_configurations)
foreach(_configuration ${_configurations})
get_target_property(_imported_location ${_library} IMPORTED_LOCATION_${_configuration})
list(APPEND _libraries ${_imported_location})
endforeach()
endif()
endif()
endforeach()
endif()

process_feature(VTK
LIBRARIES REQUIRED ${_libraries}
LIBRARIES REQUIRED _libraries
INCLUDE_DIRS REQUIRED VTK_INCLUDE_DIR
CLEAR VTK_INCLUDE_DIR VTK_LIBRARIES ${_libraries}
CLEAR VTK_INCLUDE_DIR VTK_LIBRARIES _libraries
)
81 changes: 81 additions & 0 deletions include/deal.II/vtk/utilities.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,81 @@
// ---------------------------------------------------------------------
//
// Copyright (C) 2023 by the deal.II authors
//
// This file is part of the deal.II library.
//
// The deal.II library is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE.md at
// the top level directory of deal.II.
//
// ---------------------------------------------------------------------

#ifndef dealii_vtk_utilities_h
#define dealii_vtk_utilities_h

#include <deal.II/base/config.h>

#include <deal.II/base/exceptions.h>

#ifdef DEAL_II_WITH_VTK

# include <vtkDoubleArray.h>

DEAL_II_NAMESPACE_OPEN

/**
* Interface to the Visualization Toolkit (VTK).
*
* VTK is an open-source, freely available software system for 3D computer
* graphics, modeling, image processing, volume rendering, scientific
* visualization, and 2D plotting.
* It supports a wide variety of visualization algorithms and advanced
* modeling techniques, and it takes advantage of both threaded and distributed
* memory parallel processing for speed and scalability, respectively.
*
* You can learn more about the VTK library at https://vtk.org/
*/
namespace VTKWrappers
{
/**
* Convert from a deal.II Point to a VTK double array.
*
* @tparam dim Dimension of the point
* @param [in] p An input deal.II Point<dim>
* @return A VTK smart pointer to the data array.
*/
template <int dim>
inline vtkSmartPointer<vtkDoubleArray>
dealii_point_to_vtk_array(const dealii::Point<dim> &p);

# ifndef DOXYGEN
// Template implementations

template <int dim>
inline vtkSmartPointer<vtkDoubleArray>
dealii_point_to_vtk_array(const dealii::Point<dim> &p)
{
vtkSmartPointer<vtkDoubleArray> p_vtk =
vtkSmartPointer<vtkDoubleArray>::New();

p_vtk->SetNumberOfComponents(dim);
p_vtk->SetNumberOfTuples(1);

for (int d = 0; d < dim; ++d)
p_vtk->FillComponent(d, p[d]);

return p_vtk;
}

# endif

} // namespace VTKWrappers

DEAL_II_NAMESPACE_CLOSE

#endif
#endif
6 changes: 6 additions & 0 deletions tests/vtk/CMakeLists.txt
Original file line number Diff line number Diff line change
@@ -0,0 +1,6 @@
cmake_minimum_required(VERSION 3.13.4)
include(../setup_testsubproject.cmake)
project(testsuite CXX)
if(DEAL_II_WITH_VTK)
deal_ii_pickup_tests()
endif()
39 changes: 39 additions & 0 deletions tests/vtk/vtk_01.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,39 @@
// ---------------------------------------------------------------------
//
// Copyright (C) 2023 by the deal.II authors
//
// This file is part of the deal.II library.
//
// The deal.II library is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE.md at
// the top level directory of deal.II.
//
// ---------------------------------------------------------------------

// Test conversion from deal.II Point to VTK double array

#include <deal.II/base/point.h>

#include <deal.II/vtk/utilities.h>

#include "../tests.h"

int
main()
{
initlog();

// Test conversion from deal.II Point to VTK array.
{
const Point<3> point_dealii(1.0, 2.0, 3.0);
const auto point_vtk = VTKWrappers::dealii_point_to_vtk_array(point_dealii);
deallog << "VTK Point: " << point_vtk->GetComponent(0, 0) << ", "
<< point_vtk->GetComponent(0, 1) << ", "
<< point_vtk->GetComponent(0, 2) << std::endl;
}

return 0;
}
2 changes: 2 additions & 0 deletions tests/vtk/vtk_01.output
Original file line number Diff line number Diff line change
@@ -0,0 +1,2 @@

DEAL::VTK Point: 1.00000, 2.00000, 3.00000

0 comments on commit a202b0f

Please sign in to comment.