FTK is a library that provides building blocks for feature tracking algorithms in scientific datasets. FTK is header-only and written in C++17.
- Hypermesh: data structures for high-dimensional meshes and mesh elements including n-simplices, n-cubes, and n-prisms; utilities to generalize given 2D/3D structured/unstructured meshes into 3D/4D spacetime meshes
- Numeric: root-find algorithms for inverse interpolations and parallel vector operators in n-simplices, n-cubes, and simplex-prisms; lightweight linear algebra utilities to support root-finding
- CCL: connected component labeling algorithm for building feature tracking algorithms
- Geometry: utilities to transform connect components to geometry for visualization and analysis
- Tracking graph: data structures to record births, deaths, merges, and splits of features; visualization algorithms for tracking graphs
- IO: interfaces to stage multiple timesteps of the inputs and to store outputs in key-value stores (LevelDB, RocksDB, and Mochi), file systems, and in situ staging areas
You may include FTK headers and call FTK functions directly from your C++ code, because FTK is header-only.
$ git clone https://github.com/hguo/ftk $FTK_SOURCE_DIR
You may install FTK to an designated path. The installation will also generate FTKConfig.cmake in the installation path, such that you can use find_package(FTK)
to find and use FTK in your CMakeLists.txt
$ git clone https://github.com/hguo/ftk $FTK_SOURCE_DIR
$ mkdir $FTK_SOURCE_DIR/build && cd $FTK_SOURCE_DIR/build
$ cmake .. -DCMAKE_INSTALL_PREFIX=$FTK_INSTALL_DIR
$ make install
The installed files are organized as follows:
$ tree $FTK_INSTALL_DIR
.
├── include
│ ├── ftk
│ │ ├── algorithms
│ │ │ ├── bfs.hh
...
│ └── hypermesh
│ ├── ndarray.hh
│ ├── regular_mesh.hh
...
└── lib
└── cmake
└── FTKConfig.cmake
You may use the FTK installation in your own CMakeLists.txt file:
find_package(FTK REQUIRED)
include_directories (${FTK_INCLUDE_DIR})
When you configure your build, please specify FTK_DIR with CMake:
$ cmake -DFTK_DIR=$FTK_INSTALL_DIR/lib/cmake
FTK currently provides three examples, including 2D critical point tracking, 3D critical point tracking, and a demo for the tracking graph.
FTK examples depend on the following packages
- VTK >= 8.2, optional, or required if you would like to visualize trajectories with VTK
- Qt5 >= 5.11.2, optional, or required if you would like to visualize trajectories with Qt/OpenGL windows
- NetCDF >= 4.6.1, optional, or requried if you would like to load a 3D scalar
$ cd $FTK_SOURCE_DIR/build
$ cmake .. -DFTK_BUILD_EXAMPLES=1
$ make
The followings are the guide to run the 2D critical point tracking example. We aim to track local maximum over time for the following scalar field:
where
In this example, local maximum are defined as the loci of points that
You can run the example with following command, and the executable prints the trajectory:
$ ./examples/critical_point_tracking_2d/ex_critical_point_tracking_2d
We found 16 trajectories:
--Curve 0:
---x=(2.000000, 26.972244), t=5.869000, val=0.955384
---x=(2.056124, 26.857716), t=5.857716, val=0.958166
...
--Curve 15:
---x=(125.000000, 122.414673), t=1.241837, val=0.884089
---x=(124.709686, 122.709686), t=1.187109, val=0.912740
...
If you built the example with Qt, you can use the Qt/OpenGL based window to explore the trajectory. You can rotate the space with mouse and increase/decrease timesteps with the left/right keys. The color encodes the ID of different trajectories. Please use —qt
argument to run the executable, and the example screenshot is as follows.
$ ./examples/critical_point_tracking_2d/ex_critical_point_tracking_2d --qt
If you built the example with VTK, you may also explore the trajectories with a VTK window by adding —vtk
argument:
$ ./examples/critical_point_tracking_2d/ex_critical_point_tracking_2d --vtk
- vortexfinder2: Vortex finder for time-dependent Ginzburg-Landau superconductor simulation data