![]() |
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Module providing an F90 interface to the ZFP compression library: http://computation.llnl.gov/projects/floating-point-compression In addition it provides simple threaded (de)compression routines. Important: This module uses C-like 0-based indexing!
Derived types and interfaces | |
| type | sll_t_compressed_buffer |
| data structure to support threaded ZFP compression and decompression More... | |
Functions/Subroutines | |
| subroutine | deflate_buffer_real64 (buf, comp, n_doubles, n_threads) |
| compress buffer More... | |
| subroutine | inflate_buffer_real64 (buf, comp, n_threads) |
| decompress buffer More... | |
| subroutine | deallocate_compressed_buffer_obj (comp) |
| subroutine | allocate_compressed_buffer_index_arrays (comp, n_slices) |
| subroutine | print_compression_information (comp, verbose) |
| subroutine | set_compression_precision (prec) |
| integer function | concatenate_index_arrays (comp, array) |
| allocate array, copy indices from comp into array, return More... | |
| subroutine | decatenate_index_arrays (comp, array) |
| double precision function | get_time () |
Variables | |
| integer, parameter | zfp_blocksize = 64 |
| integer | zfp_precision = 32 |
| type sll_m_compression::sll_t_compressed_buffer |
data structure to support threaded ZFP compression and decompression
Definition at line 65 of file sll_m_compression.F90.
| subroutine sll_m_compression::allocate_compressed_buffer_index_arrays | ( | type(sll_t_compressed_buffer), intent(inout) | comp, |
| integer, intent(in) | n_slices | ||
| ) |
| integer function sll_m_compression::concatenate_index_arrays | ( | type(sll_t_compressed_buffer) | comp, |
| integer, dimension(:), pointer | array | ||
| ) |
allocate array, copy indices from comp into array, return
Definition at line 296 of file sll_m_compression.F90.
| subroutine sll_m_compression::deallocate_compressed_buffer_obj | ( | type(sll_t_compressed_buffer), intent(inout) | comp | ) |
| subroutine sll_m_compression::decatenate_index_arrays | ( | type(sll_t_compressed_buffer) | comp, |
| integer, dimension(:), pointer | array | ||
| ) |
Definition at line 330 of file sll_m_compression.F90.
| subroutine sll_m_compression::deflate_buffer_real64 | ( | real(kind=f64), dimension(0:), intent(in), target | buf, |
| type(sll_t_compressed_buffer), intent(inout) | comp, | ||
| integer, intent(in), optional | n_doubles, | ||
| integer, intent(in), optional | n_threads | ||
| ) |
compress buffer
Definition at line 88 of file sll_m_compression.F90.
| double precision function sll_m_compression::get_time |
Definition at line 362 of file sll_m_compression.F90.
| subroutine sll_m_compression::inflate_buffer_real64 | ( | real(kind=f64), dimension(0:), intent(inout), target | buf, |
| type(sll_t_compressed_buffer), intent(inout) | comp, | ||
| integer, intent(in), optional | n_threads | ||
| ) |
decompress buffer
Definition at line 180 of file sll_m_compression.F90.
| subroutine sll_m_compression::print_compression_information | ( | type(sll_t_compressed_buffer), intent(in) | comp, |
| logical, intent(in), optional | verbose | ||
| ) |
| subroutine sll_m_compression::set_compression_precision | ( | integer, intent(in) | prec | ) |
| integer, parameter zfp_blocksize = 64 |
Definition at line 61 of file sll_m_compression.F90.
| integer zfp_precision = 32 |
Definition at line 62 of file sll_m_compression.F90.
1.9.1