Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Support VTK 9.x #152

Merged
merged 11 commits into from
Aug 16, 2023
2 changes: 2 additions & 0 deletions .github/workflows/portability.yml
Original file line number Diff line number Diff line change
Expand Up @@ -30,6 +30,8 @@ jobs:
extra-build-args: ""
- image: chaste/runner:jammy.portability-05
extra-build-args: ""
- image: chaste/runner:jammy.portability-06
extra-build-args: ""
- image: chaste/runner:jammy.portability-03
extra-build-args: "-DChaste_USE_CVODE=OFF -DChaste_USE_VTK=OFF"

Expand Down
65 changes: 45 additions & 20 deletions CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -194,6 +194,20 @@ if (WIN32 OR CYGWIN)
endif (WIN32 OR CYGWIN)


###############################################
#### Create interface library for dependencies
###############################################

# See https://github.com/Chaste/Chaste/issues/36 re supporting VTK 9.x.
#
# Approach is to create a interface library (with no source files) that all Chaste dependencies are linked to. This
# library is then linked against all Chaste libraries and executables.
#
# This allows modern CMake imported targets to be used for our dependencies, rather than populating lists of libaries
# and include directories as we currently to. For now, this will only affect VTK 9.x, but should be expanded to all
# dependencies in the future.
add_library(Chaste_COMMON_DEPS INTERFACE)

################################
# FIND THIRD PARTY LIBRARIES #
################################
Expand Down Expand Up @@ -322,33 +336,40 @@ if (Chaste_AUTO_INSTALL_DEPS)
endforeach ()
endif ()


################################
#### Find VTK
################################
if (Chaste_USE_VTK)
# First pass will identify the version number
find_package (VTK REQUIRED)
message (STATUS "VTK version: ${VTK_MAJOR_VERSION}.${VTK_MINOR_VERSION}")

# Depending on the VTK version, we need different components. 5.x require the following
if (VTK_MAJOR_VERSION EQUAL 5)
find_package(VTK COMPONENTS vtkIO vtkCommon vtkGraphics z REQUIRED)
# 6.0 and 6.1 require the following
elseif ((VTK_MAJOR_VERSION EQUAL 6) AND (VTK_MINOR_VERSION LESS 2))
find_package(
VTK COMPONENTS vtkCommonCore vtkCommonDataModel vtkFiltersCore vtkFiltersGeneral vtkFiltersGeneric
vtkFiltersGeometry vtkFiltersModeling vtkFiltersSources vtkIOCore vtkIOGeometry vtkIOLegacy vtkIOXML
REQUIRED
)
# 6.2 and up, up to and including 8.x, require the following
elseif (VTK_MAJOR_VERSION LESS 9)
# The complexity here is for two reasons:
# 1. VTK has renamed its components several times
# 2. From VTK 9, a simple call to find_package(VTK) will try to find QT, which is not a Chaste dependency
#
# We first find the VTK common library, which can have one of two names: this will populate the VTK version
# number. This is then used to find the specific list of required VTK components.

# List of possible base modules across different VTK versions
set(POSSIBLE_VTK_BASE_MODULES vtkCommonCore CommonCore)

foreach (BASE_MODULE ${POSSIBLE_VTK_BASE_MODULES})
find_package(VTK COMPONENTS ${BASE_MODULE} QUIET)
if (VTK_FOUND)
message (STATUS "VTK version: ${VTK_MAJOR_VERSION}.${VTK_MINOR_VERSION}")
break()
endif ()
endforeach ()

if (NOT ${VTK_FOUND})
message(FATAL_ERROR "Did not find VTK")
endif ()

# Depending on the VTK version, we need different components.
if (VTK_MAJOR_VERSION LESS 9)
find_package(
VTK COMPONENTS vtkCommonCore vtkCommonDataModel vtkFiltersCore vtkFiltersGeneral vtkFiltersGeneric
vtkFiltersGeometry vtkFiltersModeling vtkFiltersSources vtkIOCore vtkIOGeometry vtkIOLegacy
vtkIOParallelXML vtkIOXML REQUIRED
)
# VTK 9.0 and above require the following
# VTK 9.0 and above require the following
else ()
find_package(
VTK COMPONENTS CommonCore CommonDataModel FiltersCore FiltersGeneral FiltersGeneric FiltersGeometry
Expand Down Expand Up @@ -503,8 +524,12 @@ list (APPEND Chaste_LINK_LIBRARIES "${MPI_CXX_LIBRARIES}")

# make sure VTK libraries added after HDF5 so VTK's link with HDF5 isn't used
if (Chaste_USE_VTK)
list (APPEND Chaste_INCLUDES "${VTK_INCLUDE_DIRS}")
list (APPEND Chaste_LINK_LIBRARIES "${VTK_LIBRARIES}")
if (VTK_MAJOR_VERSION LESS 9)
list(APPEND Chaste_INCLUDES "${VTK_INCLUDE_DIRS}")
list(APPEND Chaste_LINK_LIBRARIES "${VTK_LIBRARIES}")
else ()
target_link_libraries(Chaste_COMMON_DEPS INTERFACE ${VTK_LIBRARIES})
endif ()
endif ()


Expand Down
4 changes: 2 additions & 2 deletions cell_based/test/population/TestVertexBasedCellPopulation.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -1607,7 +1607,7 @@ class TestVertexBasedCellPopulation : public AbstractCellBasedTestSuite
TS_ASSERT(vtk_file.Exists());

// Check that we have 144 (6*4*6) edges and 68 cells
ifstream vtk_file_stream;
std::ifstream vtk_file_stream;
vtk_file_stream.open (vtk_file.GetAbsolutePath());
std::stringstream vtk_file_string_buffer;
vtk_file_string_buffer << vtk_file_stream.rdbuf();
Expand Down Expand Up @@ -1678,7 +1678,7 @@ class TestVertexBasedCellPopulation : public AbstractCellBasedTestSuite
TS_ASSERT(vtk_file.Exists());

// Check that we have 144 (6*4*6) edges and 68 cells
ifstream vtk_file_stream;
std::ifstream vtk_file_stream;
vtk_file_stream.open (vtk_file.GetAbsolutePath());
std::stringstream vtk_file_string_buffer;
vtk_file_string_buffer << vtk_file_stream.rdbuf();
Expand Down
4 changes: 4 additions & 0 deletions cmake/Modules/ChasteMacros.cmake
Original file line number Diff line number Diff line change
Expand Up @@ -147,6 +147,7 @@ macro(Chaste_ADD_TEST _testTargetName _filename)
set_source_files_properties("${_test_real_output_filename}" PROPERTIES GENERATED true)

add_executable(${exeTargetName} "${_test_real_output_filename}" ${_filename} ${ARGN})
target_link_libraries(${exeTargetName} PRIVATE Chaste_COMMON_DEPS)
endif()

set(test_exe $<TARGET_FILE:${exeTargetName}>)
Expand Down Expand Up @@ -318,6 +319,7 @@ macro(Chaste_DO_COMMON component)
# Make component library, if component contains any source files
if (NOT Chaste_${component}_SOURCES STREQUAL "")
add_library(chaste_${component} ${Chaste_${component}_SOURCES} ${ARGN})
target_link_libraries(chaste_${component} PRIVATE Chaste_COMMON_DEPS)
if (BUILD_SHARED_LIBS)
target_link_libraries(chaste_${component} LINK_PUBLIC ${Chaste_LIBRARIES})
set(static_extension "a")
Expand Down Expand Up @@ -429,6 +431,7 @@ macro(Chaste_DO_APPS_COMMON component)
set(component_library )
endif()
add_executable(${appName} ${app})
target_link_libraries(${appName} PRIVATE Chaste_COMMON_DEPS)
#set_target_properties(${appName} PROPERTIES RUNTIME_OUTPUT_DIRECTORY ${CMAKE_CURRENT_BINARY_DIR}/src)

if (NOT ${component} STREQUAL "")
Expand Down Expand Up @@ -528,6 +531,7 @@ macro(Chaste_DO_TEST_COMMON component)
file(GLOB_RECURSE test_sources RELATIVE ${CMAKE_CURRENT_SOURCE_DIR} *.cpp)
if(test_sources)
add_library(test${component} STATIC ${test_sources})
target_link_libraries(test${component} PRIVATE Chaste_COMMON_DEPS)
set(COMPONENT_LIBRARIES ${COMPONENT_LIBRARIES} test${component})
endif()

Expand Down
1 change: 1 addition & 0 deletions global/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -57,6 +57,7 @@ int main( )
")

add_executable(timekeeper "${generate_dir}/timekeeper.cpp")
target_link_libraries(timekeeper PRIVATE Chaste_COMMON_DEPS)

###################################################################
# Setup command to generate Version and BuildInfo at build-time #
Expand Down
12 changes: 7 additions & 5 deletions mesh/src/reader/VtkMeshReader.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -103,11 +103,16 @@ void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::CommonConstructor()
mVtkCellType = VTK_LINE;
}

//Determine if we have multiple cell types - such as cable elements in addition to tets/triangles
#if (VTK_MAJOR_VERSION >= 9 && VTK_MINOR_VERSION >= 2)
const auto num_distinct_cell_types = static_cast<unsigned>(mpVtkUnstructuredGrid->GetDistinctCellTypesArray()->GetNumberOfTuples());
#else // VTK older than 9.2
vtkCellTypes* cell_types = vtkCellTypes::New();
mpVtkUnstructuredGrid->GetCellTypes(cell_types);
const auto num_distinct_cell_types = static_cast<unsigned>(cell_types->GetNumberOfTypes());
cell_types->Delete();
#endif

if (cell_types->GetNumberOfTypes() > 1)
if (num_distinct_cell_types > 1u)
{
mNumCableElementAttributes = 1;
for (unsigned cell_id = 0; cell_id < num_cells; ++cell_id)
Expand All @@ -133,9 +138,6 @@ void VtkMeshReader<ELEMENT_DIM,SPACE_DIM>::CommonConstructor()
mNumElements = num_cells;
}

cell_types->Delete();


// Extract the surface faces
if (ELEMENT_DIM == 2u)
{
Expand Down
34 changes: 21 additions & 13 deletions mesh/test/reader/TestVtkMeshReader.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -73,7 +73,18 @@ typedef VtkMeshReader<3,3> MESH_READER3;

class TestVtkMeshReader : public CxxTest::TestSuite
{
//Requires "sudo aptitude install libvtk5-dev" or similar
private:
bool IsRotation(const std::vector<unsigned>& sequence, const std::vector<unsigned>& target)
{
// Double the input sequence
std::vector<unsigned> doubled_sequence = sequence;
doubled_sequence.insert(doubled_sequence.end(), sequence.begin(), sequence.end());

// Search for the target sequence in the doubled sequence
return std::search(doubled_sequence.begin(), doubled_sequence.end(), target.begin(), target.end()) != doubled_sequence.end();
}

//Requires "sudo aptitude install libvtk5-dev" or similar
public:

/**
Expand Down Expand Up @@ -202,27 +213,24 @@ class TestVtkMeshReader : public CxxTest::TestSuite
TS_ASSERT_EQUALS(mesh_reader.GetNumFaces(), 20u);
TS_ASSERT_EQUALS(mesh_reader.GetNumEdges(), 20u);

// We don't care which order the nodes are in the face data, as long as they are correct up to a rotation.
// For instance, here, first_face_data.NodeIndices could either be {11, 3, 0}, {0, 11, 3} or {3, 0, 11}.
// VTK 9.2 started producing a different ordering, which necessitated updating this test.
// See https://github.com/Chaste/Chaste/issues/36

ElementData first_face_data = mesh_reader.GetNextFaceData();
TS_ASSERT_EQUALS(first_face_data.NodeIndices[0], 11u);
TS_ASSERT_EQUALS(first_face_data.NodeIndices[1], 3u);
TS_ASSERT_EQUALS(first_face_data.NodeIndices[2], 0u);
TS_ASSERT(IsRotation(first_face_data.NodeIndices, {11u, 3u, 0u}))

ElementData next_face_data = mesh_reader.GetNextFaceData();
TS_ASSERT_EQUALS(next_face_data.NodeIndices[0], 3u);
TS_ASSERT_EQUALS(next_face_data.NodeIndices[1], 8u);
TS_ASSERT_EQUALS(next_face_data.NodeIndices[2], 0u);
TS_ASSERT(IsRotation(next_face_data.NodeIndices, {3u, 8u, 0u}))

mesh_reader.Reset();

first_face_data = mesh_reader.GetNextEdgeData();
TS_ASSERT_EQUALS(first_face_data.NodeIndices[0], 11u);
TS_ASSERT_EQUALS(first_face_data.NodeIndices[1], 3u);
TS_ASSERT_EQUALS(first_face_data.NodeIndices[2], 0u);
TS_ASSERT(IsRotation(first_face_data.NodeIndices, {11u, 3u, 0u}))

next_face_data = mesh_reader.GetNextEdgeData();
TS_ASSERT_EQUALS(next_face_data.NodeIndices[0], 3u);
TS_ASSERT_EQUALS(next_face_data.NodeIndices[1], 8u);
TS_ASSERT_EQUALS(next_face_data.NodeIndices[2], 0u);
TS_ASSERT(IsRotation(next_face_data.NodeIndices, {3u, 8u, 0u}))

for (unsigned face=2; face<mesh_reader.GetNumFaces(); face++)
{
Expand Down