Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
9 changes: 8 additions & 1 deletion httomolibgpu/__init__.py
Original file line number Diff line number Diff line change
Expand Up @@ -11,5 +11,12 @@
remove_all_stripe,
)

from httomolibgpu.recon.algorithm import FBP2d_astra, FBP, LPRec, SIRT, CGLS
from httomolibgpu.recon.algorithm import (
FBP2d_astra,
FBP3d_tomobar,
LPRec3d_tomobar,
SIRT3d_tomobar,
CGLS3d_tomobar,
)

from httomolibgpu.recon.rotation import find_center_vo, find_center_360, find_center_pc
24 changes: 12 additions & 12 deletions httomolibgpu/recon/algorithm.py
Original file line number Diff line number Diff line change
Expand Up @@ -43,10 +43,10 @@

__all__ = [
"FBP2d_astra",
"FBP",
"SIRT",
"CGLS",
"LPRec",
"FBP3d_tomobar",
"LPRec3d_tomobar",
"SIRT3d_tomobar",
"CGLS3d_tomobar",
]

input_data_axis_labels = ["angles", "detY", "detX"] # set the labels of the input data
Expand Down Expand Up @@ -130,13 +130,13 @@ def FBP2d_astra(


## %%%%%%%%%%%%%%%%%%%%%%% FBP reconstruction %%%%%%%%%%%%%%%%%%%%%%%%%%%% ##
def FBP(
def FBP3d_tomobar(
data: cp.ndarray,
angles: np.ndarray,
center: Optional[float] = None,
filter_freq_cutoff: float = 0.35,
recon_size: Optional[int] = None,
recon_mask_radius: float = 0.95,
recon_mask_radius: Optional[float] = 0.95,
neglog: bool = False,
gpu_id: int = 0,
) -> cp.ndarray:
Expand All @@ -158,7 +158,7 @@ def FBP(
recon_size : int, optional
The [recon_size, recon_size] shape of the reconstructed slice in pixels.
By default (None), the reconstructed size will be the dimension of the horizontal detector.
recon_mask_radius: float
recon_mask_radius: float, optional
The radius of the circular mask that applies to the reconstructed slice in order to crop
out some undesirable artifacts. The values outside the given diameter will be set to zero.
It is recommended to keep the value in the range [0.7-1.0].
Expand All @@ -171,7 +171,7 @@ def FBP(
Returns
-------
cp.ndarray
The FBP reconstructed volume as a CuPy array.
FBP reconstructed volume as a CuPy array.
"""
RecToolsCP = _instantiate_direct_recon_class(
data, angles, center, recon_size, gpu_id
Expand All @@ -188,7 +188,7 @@ def FBP(


## %%%%%%%%%%%%%%%%%%%%%%% LPRec %%%%%%%%%%%%%%%%%%%%%%%%%%%% ##
def LPRec(
def LPRec3d_tomobar(
data: cp.ndarray,
angles: np.ndarray,
center: Optional[float] = None,
Expand All @@ -212,7 +212,7 @@ def LPRec(
recon_size : int, optional
The [recon_size, recon_size] shape of the reconstructed slice in pixels.
By default (None), the reconstructed size will be the dimension of the horizontal detector.
recon_mask_radius: float
recon_mask_radius: float, optional
The radius of the circular mask that applies to the reconstructed slice in order to crop
out some undesirable artifacts. The values outside the given diameter will be set to zero.
It is recommended to keep the value in the range [0.7-1.0].
Expand All @@ -237,7 +237,7 @@ def LPRec(


## %%%%%%%%%%%%%%%%%%%%%%% SIRT reconstruction %%%%%%%%%%%%%%%%%%%%%%%%%%%% ##
def SIRT(
def SIRT3d_tomobar(
data: cp.ndarray,
angles: np.ndarray,
center: Optional[float] = None,
Expand Down Expand Up @@ -302,7 +302,7 @@ def SIRT(


## %%%%%%%%%%%%%%%%%%%%%%% CGLS reconstruction %%%%%%%%%%%%%%%%%%%%%%%%%%%% ##
def CGLS(
def CGLS3d_tomobar(
data: cp.ndarray,
angles: np.ndarray,
center: Optional[float] = None,
Expand Down
43 changes: 21 additions & 22 deletions tests/test_recon/test_algorithm.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,10 +3,10 @@
from httomolibgpu.prep.normalize import normalize as normalize_cupy
from httomolibgpu.recon.algorithm import (
FBP2d_astra,
FBP,
LPRec,
SIRT,
CGLS,
FBP3d_tomobar,
LPRec3d_tomobar,
SIRT3d_tomobar,
CGLS3d_tomobar,
)
from numpy.testing import assert_allclose
import time
Expand Down Expand Up @@ -34,9 +34,8 @@ def test_reconstruct_FBP_2d_astra(data, flats, darks, ensure_clean_memory):
assert recon_data.dtype == np.float32
assert recon_data.shape == (recon_size, 128, recon_size)


def test_reconstruct_FBP_1(data, flats, darks, ensure_clean_memory):
recon_data = FBP(
def test_reconstruct_FBP3d_tomobar_1(data, flats, darks, ensure_clean_memory):
recon_data = FBP3d_tomobar(
normalize_cupy(data, flats, darks, cutoff=10, minus_log=True),
np.linspace(0.0 * np.pi / 180.0, 180.0 * np.pi / 180.0, data.shape[0]),
79.5,
Expand All @@ -53,8 +52,8 @@ def test_reconstruct_FBP_1(data, flats, darks, ensure_clean_memory):
assert recon_data.shape == (160, 128, 160)


def test_reconstruct_FBP_1_neglog(data, flats, darks, ensure_clean_memory):
recon_data = FBP(
def test_reconstruct_FBP3d_tomobar_1_neglog(data, flats, darks, ensure_clean_memory):
recon_data = FBP3d_tomobar(
normalize_cupy(data, flats, darks, cutoff=10, minus_log=False),
np.linspace(0.0 * np.pi / 180.0, 180.0 * np.pi / 180.0, data.shape[0]),
79.5,
Expand All @@ -72,8 +71,8 @@ def test_reconstruct_FBP_1_neglog(data, flats, darks, ensure_clean_memory):
assert recon_data.shape == (160, 128, 160)


def test_reconstruct_FBP_2(data, flats, darks, ensure_clean_memory):
recon_data = FBP(
def test_reconstruct_FBP3d_tomobar_2(data, flats, darks, ensure_clean_memory):
recon_data = FBP3d_tomobar(
normalize_cupy(data, flats, darks, cutoff=20.5, minus_log=False),
np.linspace(5.0 * np.pi / 360.0, 180.0 * np.pi / 360.0, data.shape[0]),
15.5,
Expand All @@ -90,8 +89,8 @@ def test_reconstruct_FBP_2(data, flats, darks, ensure_clean_memory):
assert recon_data.dtype == np.float32


def test_reconstruct_FBP_3(data, flats, darks, ensure_clean_memory):
recon_data = FBP(
def test_reconstruct_FBP3d_tomobar_3(data, flats, darks, ensure_clean_memory):
recon_data = FBP3d_tomobar(
normalize_cupy(data, flats, darks, cutoff=20.5, minus_log=False),
np.linspace(5.0 * np.pi / 360.0, 180.0 * np.pi / 360.0, data.shape[0]),
79, # center
Expand All @@ -109,8 +108,8 @@ def test_reconstruct_FBP_3(data, flats, darks, ensure_clean_memory):
assert recon_data.shape == (210, 128, 210)


def test_reconstruct_LPREC_1(data, flats, darks, ensure_clean_memory):
recon_data = LPRec(
def test_reconstruct_LPRec3d_tomobar_1(data, flats, darks, ensure_clean_memory):
recon_data = LPRec3d_tomobar(
data=normalize_cupy(data, flats, darks, cutoff=10, minus_log=True),
angles=np.linspace(0.0 * np.pi / 180.0, 180.0 * np.pi / 180.0, data.shape[0]),
center=79.5,
Expand All @@ -124,9 +123,9 @@ def test_reconstruct_LPREC_1(data, flats, darks, ensure_clean_memory):
assert recon_data.shape == (130, 128, 130)


def test_reconstruct_SIRT(data, flats, darks, ensure_clean_memory):
def test_reconstruct_SIRT3d_tomobar(data, flats, darks, ensure_clean_memory):
objrecon_size = data.shape[2]
recon_data = SIRT(
recon_data = SIRT3d_tomobar(
normalize_cupy(data, flats, darks, cutoff=10, minus_log=True),
np.linspace(0.0 * np.pi / 180.0, 180.0 * np.pi / 180.0, data.shape[0]),
79.5,
Expand All @@ -141,9 +140,9 @@ def test_reconstruct_SIRT(data, flats, darks, ensure_clean_memory):
assert recon_data.dtype == np.float32


def test_reconstruct_CGLS(data, flats, darks, ensure_clean_memory):
def test_reconstruct_CGLS3d_tomobar(data, flats, darks, ensure_clean_memory):
objrecon_size = data.shape[2]
recon_data = CGLS(
recon_data = CGLS3d_tomobar(
normalize_cupy(data, flats, darks, cutoff=10, minus_log=True),
np.linspace(0.0 * np.pi / 180.0, 180.0 * np.pi / 180.0, data.shape[0]),
79.5,
Expand All @@ -159,7 +158,7 @@ def test_reconstruct_CGLS(data, flats, darks, ensure_clean_memory):


@pytest.mark.perf
def test_FBP_performance(ensure_clean_memory):
def test_FBP3d_tomobar_performance(ensure_clean_memory):
dev = cp.cuda.Device()
data_host = np.random.random_sample(size=(1801, 5, 2560)).astype(np.float32) * 2.0
data = cp.asarray(data_host, dtype=np.float32)
Expand All @@ -168,13 +167,13 @@ def test_FBP_performance(ensure_clean_memory):
filter_freq_cutoff = 1.1

# cold run first
FBP(data, angles, cor, filter_freq_cutoff)
FBP3d_tomobar(data, angles, cor, filter_freq_cutoff)
dev.synchronize()

start = time.perf_counter_ns()
nvtx.RangePush("Core")
for _ in range(10):
FBP(data, angles, cor)
FBP3d_tomobar(data, angles, cor)
nvtx.RangePop()
dev.synchronize()
duration_ms = float(time.perf_counter_ns() - start) * 1e-6 / 10
Expand Down
35 changes: 17 additions & 18 deletions zenodo-tests/test_recon/test_algorithm.py
Original file line number Diff line number Diff line change
Expand Up @@ -8,8 +8,8 @@
from httomolibgpu.prep.normalize import normalize
from httomolibgpu.recon.algorithm import (
FBP2d_astra,
FBP,
LPRec,
FBP3d_tomobar,
LPRec3d_tomobar,
)
from httomolibgpu.misc.morph import sino_360_to_180
from numpy.testing import assert_allclose
Expand All @@ -18,7 +18,6 @@
from cupy.cuda import nvtx
from conftest import force_clean_gpu_memory


def test_reconstruct_FBP2d_astra_i12_dataset1(i12_dataset1):
force_clean_gpu_memory()
projdata = i12_dataset1[0]
Expand All @@ -44,9 +43,9 @@ def test_reconstruct_FBP2d_astra_i12_dataset1(i12_dataset1):
assert_allclose(np.sum(recon_data), 84673.68, rtol=1e-07, atol=1e-6)
assert recon_data.dtype == np.float32
assert recon_data.shape == (2560, 50, 2560)



def test_reconstruct_FBP_i12_dataset1(i12_dataset1):
def test_reconstruct_FBP3d_tomobar_i12_dataset1(i12_dataset1):
force_clean_gpu_memory()
projdata = i12_dataset1[0]
angles = i12_dataset1[1]
Expand All @@ -58,7 +57,7 @@ def test_reconstruct_FBP_i12_dataset1(i12_dataset1):
del flats, darks, projdata
force_clean_gpu_memory()

recon_data = FBP(
recon_data = FBP3d_tomobar(
data_normalised,
np.deg2rad(angles),
center=1253.75,
Expand Down Expand Up @@ -92,13 +91,13 @@ def test_reconstruct_LP_REC_i13_dataset1(i13_dataset1):
# GPU archetictures older than 5.3 wont accept the data larger than
# (4096, 4096, 4096), while the newer ones can accept (16384 x 16384 x 16384)

# recon_data = FBP(
# recon_data = FBP3d_tomobar(
# stiched_data_180degrees,
# np.deg2rad(angles[0:3000]),
# center=2322,
# filter_freq_cutoff=0.35,
# )
recon_data = LPRec(
recon_data = LPRec3d_tomobar(
data=stiched_data_180degrees,
angles=np.deg2rad(angles[0:3000]),
center=2322.08,
Expand All @@ -112,7 +111,7 @@ def test_reconstruct_LP_REC_i13_dataset1(i13_dataset1):


@pytest.mark.perf
def test_FBP_performance_i13_dataset2(i13_dataset2):
def test_FBP3d_tomobar_performance_i13_dataset2(i13_dataset2):
force_clean_gpu_memory()
dev = cp.cuda.Device()
projdata = i13_dataset2[0]
Expand All @@ -126,7 +125,7 @@ def test_FBP_performance_i13_dataset2(i13_dataset2):
force_clean_gpu_memory()

# cold run first
FBP(
FBP3d_tomobar(
data_normalised,
np.deg2rad(angles),
center=1253.75,
Expand All @@ -137,7 +136,7 @@ def test_FBP_performance_i13_dataset2(i13_dataset2):
start = time.perf_counter_ns()
nvtx.RangePush("Core")
for _ in range(10):
FBP(
FBP3d_tomobar(
data_normalised,
np.deg2rad(angles),
center=1286.25,
Expand All @@ -150,7 +149,7 @@ def test_FBP_performance_i13_dataset2(i13_dataset2):
assert "performance in ms" == duration_ms


def test_reconstruct_LPREC_i13_dataset2(i13_dataset2):
def test_reconstruct_LPRec3d_tomobar_i13_dataset2(i13_dataset2):
force_clean_gpu_memory()
projdata = i13_dataset2[0]
angles = i13_dataset2[1]
Expand All @@ -162,7 +161,7 @@ def test_reconstruct_LPREC_i13_dataset2(i13_dataset2):
del flats, darks, projdata
force_clean_gpu_memory()

recon_data = LPRec(
recon_data = LPRec3d_tomobar(
data=data_normalised,
angles=np.deg2rad(angles),
center=1286.25,
Expand All @@ -176,7 +175,7 @@ def test_reconstruct_LPREC_i13_dataset2(i13_dataset2):


@pytest.mark.perf
def test_LPREC_performance_i13_dataset2(i13_dataset2):
def test_LPRec3d_tomobar_performance_i13_dataset2(i13_dataset2):
dev = cp.cuda.Device()
projdata = i13_dataset2[0]
angles = i13_dataset2[1]
Expand All @@ -189,7 +188,7 @@ def test_LPREC_performance_i13_dataset2(i13_dataset2):
force_clean_gpu_memory()

# cold run first
LPRec(
LPRec3d_tomobar(
data=data_normalised,
angles=np.deg2rad(angles),
center=1286.25,
Expand All @@ -199,7 +198,7 @@ def test_LPREC_performance_i13_dataset2(i13_dataset2):
start = time.perf_counter_ns()
nvtx.RangePush("Core")
for _ in range(10):
LPRec(
LPRec3d_tomobar(
data=data_normalised,
angles=np.deg2rad(angles),
center=1286.25,
Expand All @@ -211,7 +210,7 @@ def test_LPREC_performance_i13_dataset2(i13_dataset2):
assert "performance in ms" == duration_ms


def test_reconstruct_FBP_i13_dataset3(i13_dataset3):
def test_reconstruct_FBP3d_tomobar_i13_dataset3(i13_dataset3):
force_clean_gpu_memory()
projdata = i13_dataset3[0]
angles = i13_dataset3[1]
Expand All @@ -231,7 +230,7 @@ def test_reconstruct_FBP_i13_dataset3(i13_dataset3):
# GPU archetictures older than 5.3 wont accept the data larger than
# (4096, 4096, 4096), while the newer ones can accept (16384 x 16384 x 16384)

recon_data = FBP(
recon_data = FBP3d_tomobar(
stiched_data_180degrees,
np.deg2rad(angles[0:3000]),
center=2341,
Expand Down