Skip to content
GitLab
Projects
Groups
Snippets
/
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
Menu
Open sidebar
tomotools
tomoscan
Commits
2bc2c976
Commit
2bc2c976
authored
Mar 09, 2020
by
payno
Browse files
Merge branch 'master' of gitlab.esrf.fr:tomotools/tomoscan
parents
ec0f7456
ff26348e
Pipeline
#22669
passed with stages
in 2 minutes and 12 seconds
Changes
9
Pipelines
1
Expand all
Hide whitespace changes
Inline
Side-by-side
doc/index.rst
View file @
2bc2c976
tomoscan
========
tomoscan aims to provide an unified interface to read tomography data from .edf or .hdf5 (NXtomo) acquisitions.
.. toctree::
:maxdepth: 4
...
...
tomoscan/esrf/edfscan.py
View file @
2bc2c976
...
...
@@ -22,6 +22,8 @@
#
#############################################################################*/
"""contains EDFTomoScan, class to be used with EDF acquisition"""
__authors__
=
[
"H.Payno"
]
__license__
=
"MIT"
...
...
@@ -38,7 +40,7 @@ import io
from
typing
import
Union
from
..scanbase
import
TomoScanBase
from
.utils
import
get_parameters_frm_par_or_info
,
extract_urls_from_edf
from
..unitsystem
import
m
etric
s
ystem
from
..unitsystem
.metricsystem
import
M
etric
S
ystem
from
..utils
import
docstring
import
logging
...
...
@@ -84,8 +86,25 @@ class EDFTomoScan(TomoScanBase):
self
.
__ref_on
=
None
self
.
__scan_range
=
None
self
.
_edf_n_frames
=
n_frames
self
.
__distance
=
None
self
.
__energy
=
None
self
.
update
()
@
docstring
(
TomoScanBase
.
clear_caches
)
def
clear_caches
(
self
):
self
.
_darks
=
None
self
.
_flats
=
None
self
.
_projections
=
None
self
.
__tomo_n
=
None
self
.
__ref_n
=
None
self
.
__dark_n
=
None
self
.
__dim1
=
None
self
.
__dim2
=
None
self
.
__pixel_size
=
None
self
.
__ref_on
=
None
self
.
__scan_range
=
None
@
docstring
(
TomoScanBase
.
tomo_n
)
@
property
def
tomo_n
(
self
)
->
Union
[
None
,
int
]:
...
...
@@ -116,7 +135,7 @@ class EDFTomoScan(TomoScanBase):
:rtype: float
"""
if
self
.
__pixel_size
is
None
:
self
.
__pixel_size
=
EDFTomoScan
.
get_pixel_size
(
scan
=
self
.
path
)
self
.
__pixel_size
=
EDFTomoScan
.
_
get_pixel_size
(
scan
=
self
.
path
)
return
self
.
__pixel_size
@
property
...
...
@@ -190,7 +209,7 @@ class EDFTomoScan(TomoScanBase):
def
is_abort
(
self
,
**
kwargs
)
->
bool
:
abort_file
=
os
.
path
.
basename
(
self
.
path
)
+
self
.
ABORT_FILE
abort_file
=
os
.
path
.
join
(
self
.
path
,
abort_file
)
if
'src_pattern'
in
kwargs
and
kwargs
[
'src_pattern'
is
not
None
]
:
if
'src_pattern'
in
kwargs
and
kwargs
[
'src_pattern'
]
is
not
None
:
assert
'dest_pattern'
in
kwargs
abort_file
=
abort_file
.
replace
(
kwargs
[
'src_pattern'
],
kwargs
[
'dest_pattern'
])
...
...
@@ -218,9 +237,10 @@ class EDFTomoScan(TomoScanBase):
@
docstring
(
TomoScanBase
.
update
)
def
update
(
self
):
self
.
projections
=
EDFTomoScan
.
get_proj_urls
(
self
.
path
,
n_frames
=
self
.
_edf_n_frames
)
self
.
_darks
=
EDFTomoScan
.
get_darks_url
(
self
.
path
)
self
.
_flats
=
EDFTomoScan
.
get_refs_url
(
self
.
path
)
if
self
.
path
is
not
None
:
self
.
projections
=
EDFTomoScan
.
get_proj_urls
(
self
.
path
,
n_frames
=
self
.
_edf_n_frames
)
self
.
_darks
=
EDFTomoScan
.
get_darks_url
(
self
.
path
)
self
.
_flats
=
EDFTomoScan
.
get_refs_url
(
self
.
path
)
@
docstring
(
TomoScanBase
.
load_from_dict
)
def
load_from_dict
(
self
,
desc
:
Union
[
dict
,
io
.
TextIOWrapper
]):
...
...
@@ -228,11 +248,11 @@ class EDFTomoScan(TomoScanBase):
data
=
json
.
load
(
desc
)
else
:
data
=
desc
if
not
(
self
.
_
DICT_TYPE_KEY
in
data
and
data
[
self
.
_
DICT_TYPE_KEY
]
==
self
.
_TYPE
):
if
not
(
self
.
DICT_TYPE_KEY
in
data
and
data
[
self
.
DICT_TYPE_KEY
]
==
self
.
_TYPE
):
raise
ValueError
(
'Description is not an EDFScan json description'
)
assert
self
.
_
DICT_PATH_KEY
in
data
self
.
path
=
data
[
self
.
_
DICT_PATH_KEY
]
assert
self
.
DICT_PATH_KEY
in
data
self
.
path
=
data
[
self
.
DICT_PATH_KEY
]
return
self
@
staticmethod
...
...
@@ -303,8 +323,7 @@ class EDFTomoScan(TomoScanBase):
def
extract_index
(
my_str
,
type_
):
res
=
[]
modified_str
=
copy
.
copy
(
my_str
)
while
modified_str
[
-
1
].
isdigit
():
while
modified_str
!=
""
and
modified_str
[
-
1
].
isdigit
():
res
.
append
(
modified_str
[
-
1
])
modified_str
=
modified_str
[:
-
1
]
if
len
(
res
)
==
0
:
...
...
@@ -315,7 +334,7 @@ class EDFTomoScan(TomoScanBase):
return
int
(
''
.
join
(
orignalOrder
)),
modified_str
else
:
return
float
(
'.'
.
join
((
'0'
,
''
.
join
(
orignalOrder
)))),
modified_str
_file
=
os
.
path
.
basename
(
_file
)
if
_file
.
endswith
(
'.edf'
):
name
=
_file
.
replace
(
basename
,
''
,
1
)
name
=
name
.
rstrip
(
'.edf'
)
...
...
@@ -330,7 +349,10 @@ class EDFTomoScan(TomoScanBase):
if
part_1
is
None
:
return
None
if
part_2
is
None
:
return
int
(
part_1
)
if
part_1
is
None
:
return
None
else
:
return
int
(
part_1
)
else
:
return
float
(
part_1
)
+
part_2
else
:
...
...
@@ -392,8 +414,34 @@ class EDFTomoScan(TomoScanBase):
return
d1
,
d2
@
property
@
docstring
(
TomoScanBase
.
distance
)
def
distance
(
self
)
->
Union
[
None
,
float
]:
if
self
.
__distance
is
None
:
self
.
__distance
=
EDFTomoScan
.
retrieve_information
(
self
.
path
,
None
,
"Distance"
,
type_
=
float
,
key_aliases
=
(
'distance'
,
)
)
if
self
.
__distance
is
None
:
return
None
else
:
return
self
.
__distance
*
MetricSystem
.
MILLIMETER
.
value
@
property
@
docstring
(
TomoScanBase
.
energy
)
def
energy
(
self
):
if
self
.
__energy
is
None
:
self
.
__energy
=
EDFTomoScan
.
retrieve_information
(
self
.
path
,
None
,
"Energy"
,
type_
=
float
,
key_aliases
=
(
'energy'
,
)
)
return
self
.
__energy
@
staticmethod
def
get_pixel_size
(
scan
:
str
)
->
Union
[
None
,
in
t
]:
def
_
get_pixel_size
(
scan
:
str
)
->
Union
[
None
,
floa
t
]:
if
os
.
path
.
isdir
(
scan
)
is
False
:
return
None
value
=
EDFTomoScan
.
retrieve_information
(
scan
=
scan
,
...
...
@@ -413,7 +461,7 @@ class EDFTomoScan(TomoScanBase):
# for now pixel size are stored in microns.
# We want to return them in meter
if
value
is
not
None
:
return
value
*
m
etric
s
ystem
.
micrometer
return
value
*
M
etric
S
ystem
.
MICROMETER
.
value
else
:
return
None
...
...
tomoscan/esrf/hdf5scan.py
View file @
2bc2c976
This diff is collapsed.
Click to expand it.
tomoscan/esrf/test/test_edfscan.py
View file @
2bc2c976
...
...
@@ -259,7 +259,9 @@ class TestProjections(unittest.TestCase):
self
.
assertEqual
(
len
(
scan
.
projections
),
3
)
scan
.
update
()
self
.
assertEqual
(
len
(
scan
.
projections
),
4
)
self
.
assertTrue
(
isinstance
(
scan
.
projections
[
0
],
silx
.
io
.
url
.
DataUrl
))
index_0
=
list
(
scan
.
projections
.
keys
())[
0
]
self
.
assertTrue
(
isinstance
(
scan
.
projections
[
index_0
],
silx
.
io
.
url
.
DataUrl
))
def
testProjectionWithExtraRadio
(
self
):
mock
=
MockEDF
(
scan_path
=
self
.
folder
,
n_radio
=
11
,
n_extra_radio
=
2
,
...
...
@@ -305,10 +307,6 @@ class TestScanValidatorFindFiles(unittest.TestCase):
if
os
.
path
.
isdir
(
self
.
path
):
shutil
.
rmtree
(
self
.
path
)
def
testGetRadioPaths
(
self
):
nFound
=
len
(
EDFTomoScan
.
get_proj_urls
(
self
.
path
))
self
.
assertTrue
(
nFound
==
self
.
N_RADIO
)
class
TestRadioPath
(
unittest
.
TestCase
):
"""Test static method getRadioPaths for EDFTomoScan"""
...
...
tomoscan/esrf/test/test_hdf5scan.py
View file @
2bc2c976
...
...
@@ -33,7 +33,9 @@ import shutil
import
os
import
tempfile
from
tomoscan.test.utils
import
UtilsTest
from
tomoscan.esrf.hdf5scan
import
HDF5TomoScan
from
tomoscan.esrf.hdf5scan
import
HDF5TomoScan
,
ImageKey
,
Frame
from
...unitsystem
import
metricsystem
from
silx.io.utils
import
get_data
import
numpy
...
...
@@ -56,69 +58,148 @@ class TestHDF5Scan(HDF5TestBaseClass):
"""Basic test for the hdf5 scan"""
def
setUp
(
self
)
->
None
:
super
(
TestHDF5Scan
,
self
).
setUp
()
self
.
dataset_file
=
self
.
get_dataset
(
'
data_test2.h5
'
)
self
.
scan
=
HDF5TomoScan
(
self
.
dataset_file
)
self
.
dataset_file
=
self
.
get_dataset
(
'
frm_edftomomill_twoentries.nx
'
)
self
.
scan
=
HDF5TomoScan
(
scan
=
self
.
dataset_file
)
def
testGeneral
(
self
):
"""some general on the HDF5Scan """
self
.
assertEqual
(
self
.
scan
.
master_file
,
self
.
dataset_file
)
self
.
assertEqual
(
self
.
scan
.
path
,
os
.
path
.
dirname
(
self
.
dataset_file
))
self
.
assertEqual
(
self
.
scan
.
type
,
'hdf5'
)
self
.
assertEqual
(
self
.
scan
.
entry
,
'entry0000'
)
self
.
assertEqual
(
len
(
self
.
scan
.
flats
),
42
)
self
.
assertEqual
(
len
(
self
.
scan
.
darks
),
1
)
self
.
assertEqual
(
len
(
self
.
scan
.
return_projs
),
3
)
proj_angles
=
self
.
scan
.
get_proj_angle_url
()
self
.
assertEqual
(
len
(
proj_angles
),
1500
+
3
)
self
.
assertTrue
(
90
in
proj_angles
)
self
.
assertTrue
(
24.0
in
proj_angles
)
self
.
assertTrue
(
'90.0(1)'
in
proj_angles
)
self
.
assertTrue
(
'180.0(1)'
in
proj_angles
)
self
.
assertTrue
(
179.88
not
in
proj_angles
)
url_1
=
proj_angles
[
0
]
self
.
assertTrue
(
url_1
.
is_valid
())
self
.
assertTrue
(
url_1
.
is_absolute
())
self
.
assertEquals
(
url_1
.
scheme
(),
'silx'
)
# check conversion to dict
_dict
=
self
.
scan
.
to_dict
()
scan2
=
HDF5TomoScan
.
from_dict
(
_dict
)
self
.
assertEqual
(
scan2
.
path
,
self
.
scan
.
path
)
radios_urls
=
self
.
scan
.
get_proj_angle_url
()
self
.
assertEqual
(
len
(
radios_urls
),
100
)
url_1
=
radios_urls
[
0
]
self
.
assertTrue
(
url_1
.
is_valid
())
self
.
assertFalse
(
url_1
.
is_absolute
())
self
.
assertEquals
(
url_1
.
scheme
(),
'silx'
)
self
.
assertEqual
(
scan2
.
entry
,
self
.
scan
.
entry
)
def
testFrames
(
self
):
"""Check the `frames` property which is massively used under the
HDF5TomoScan class"""
frames
=
self
.
scan
.
frames
# check some projections
proj_2
=
frames
[
24
]
self
.
assertTrue
(
isinstance
(
proj_2
,
Frame
))
self
.
assertEqual
(
proj_2
.
index
,
24
)
numpy
.
isclose
(
proj_2
.
rotation_angle
,
0.24
)
self
.
assertFalse
(
proj_2
.
is_control
)
self
.
assertEqual
(
proj_2
.
url
.
file_path
(),
self
.
scan
.
master_file
)
self
.
assertEqual
(
proj_2
.
url
.
data_path
(),
'entry0000/instrument/detector/data'
)
self
.
assertEqual
(
proj_2
.
url
.
data_slice
(),
24
)
self
.
assertEqual
(
proj_2
.
image_key
,
ImageKey
.
PROJECTION
)
self
.
assertEqual
(
get_data
(
proj_2
.
url
).
shape
,
(
20
,
20
))
# check last two non-return projection
for
frame_index
in
(
1520
,
1542
):
with
self
.
subTest
(
frame_index
=
frame_index
):
frame
=
frames
[
frame_index
]
self
.
assertTrue
(
frame
.
image_key
,
ImageKey
.
PROJECTION
)
self
.
assertFalse
(
frame
.
is_control
)
# check some darks
dark_0
=
frames
[
0
]
self
.
assertEqual
(
dark_0
.
index
,
0
)
numpy
.
isclose
(
dark_0
.
rotation_angle
,
0.0
)
self
.
assertFalse
(
dark_0
.
is_control
)
self
.
assertEqual
(
dark_0
.
url
.
file_path
(),
self
.
scan
.
master_file
)
self
.
assertEqual
(
dark_0
.
url
.
data_path
(),
'entry0000/instrument/detector/data'
)
self
.
assertEqual
(
dark_0
.
url
.
data_slice
(),
0
)
self
.
assertEqual
(
dark_0
.
image_key
,
ImageKey
.
DARK_FIELD
)
self
.
assertEqual
(
get_data
(
dark_0
.
url
).
shape
,
(
20
,
20
))
# check some refs
ref_1
=
frames
[
2
]
self
.
assertEqual
(
ref_1
.
index
,
2
)
numpy
.
isclose
(
ref_1
.
rotation_angle
,
0.0
)
self
.
assertFalse
(
ref_1
.
is_control
)
self
.
assertEqual
(
ref_1
.
url
.
file_path
(),
self
.
scan
.
master_file
)
self
.
assertEqual
(
ref_1
.
url
.
data_path
(),
'entry0000/instrument/detector/data'
)
self
.
assertEqual
(
ref_1
.
url
.
data_slice
(),
2
)
self
.
assertEqual
(
ref_1
.
image_key
,
ImageKey
.
FLAT_FIELD
)
self
.
assertEqual
(
get_data
(
ref_1
.
url
).
shape
,
(
20
,
20
))
# check some return projections
r_proj_0
=
frames
[
1543
]
self
.
assertTrue
(
isinstance
(
r_proj_0
,
Frame
))
self
.
assertEqual
(
r_proj_0
.
index
,
1543
)
numpy
.
isclose
(
r_proj_0
.
rotation_angle
,
180
)
self
.
assertTrue
(
r_proj_0
.
is_control
)
self
.
assertEqual
(
r_proj_0
.
url
.
file_path
(),
self
.
scan
.
master_file
)
self
.
assertEqual
(
r_proj_0
.
url
.
data_path
(),
'entry0000/instrument/detector/data'
)
self
.
assertEqual
(
r_proj_0
.
url
.
data_slice
(),
1543
)
self
.
assertEqual
(
r_proj_0
.
image_key
,
ImageKey
.
PROJECTION
)
self
.
assertEqual
(
get_data
(
r_proj_0
.
url
).
shape
,
(
20
,
20
))
def
testProjections
(
self
):
"""Make sure projections are valid"""
projections
=
self
.
scan
.
projections
self
.
assertEqual
(
len
(
projections
),
100
)
proj_1
=
projections
[
0
]
self
.
assertEqual
(
proj_1
.
file_path
(),
'../../../../../../users/opid19/W:/clemence/visualtomo/data_test2/tomo0001/tomo_0000.h5'
)
self
.
assertEqual
(
proj_1
.
data_slice
(),
(
'0'
,))
self
.
assertTrue
(
100
not
in
projections
)
@
unittest
.
skip
(
'no valid hdf5 acquisition defined yet'
)
self
.
assertEqual
(
len
(
self
.
scan
.
projections
),
1500
)
url_0
=
projections
[
list
(
projections
.
keys
())[
0
]]
self
.
assertEqual
(
url_0
.
file_path
(),
os
.
path
.
join
(
self
.
scan
.
master_file
))
self
.
assertEqual
(
url_0
.
data_slice
(),
22
)
def
testDark
(
self
):
n_dark
=
20
"""Make sure darks are valid"""
n_dark
=
1
self
.
assertEqual
(
self
.
scan
.
dark_n
,
n_dark
)
with
self
.
assertRaises
(
NotImplementedError
):
self
.
scan
.
darks
darks
=
self
.
scan
.
darks
self
.
assertEqual
(
len
(
darks
),
1
)
# TODO check accumulation time
def
testFlats
(
self
):
"""Make sure flats are valid"""
n_flats
=
42
flats
=
self
.
scan
.
flats
self
.
assertEqual
(
len
(
flats
),
n_flats
)
self
.
assertEqual
(
self
.
scan
.
ref_n
,
n_flats
)
with
self
.
assertRaises
(
NotImplementedError
):
flats
=
self
.
scan
.
flats
n_ref
=
21
self
.
assertEqual
(
self
.
scan
.
ref_n
,
n_ref
)
self
.
scan
.
ff_interval
# for now not implemented
# not implemented at the moment
data
=
numpy
.
arange
(
2048
*
2048
)
data
.
reshape
(
2048
,
2048
)
# not implemented at the moment
self
.
assertEqual
(
self
.
scan
.
ff_interval
,
21
)
def
testDims
(
self
):
self
.
assertEqual
(
self
.
scan
.
dim_1
,
20
)
self
.
assertEqual
(
self
.
scan
.
dim_2
,
20
)
def
testAxisUtils
(
self
):
self
.
assertEqual
(
self
.
scan
.
scan_range
,
360
)
self
.
assertEqual
(
self
.
scan
.
tomo_n
,
100
)
# self.assertEqual(self.scan.dim_1, 2048)
# self.assertEqual(self.scan.dim_2, 2048)
self
.
assertEqual
(
self
.
scan
.
scan_range
,
180
)
self
.
assertEqual
(
self
.
scan
.
tomo_n
,
1500
)
proj0_file_path
=
'../../../../../../users/opid19/W:/clemence/visualtomo/data_test2/tomo0001/tomo_0000.h5'
radios_urls_evolution
=
self
.
scan
.
get_proj_angle_url
()
self
.
assertEquals
(
len
(
radios_urls_evolution
),
1
00
)
self
.
assertEquals
(
radios_urls_evolution
[
0
].
file_path
(),
proj0_file_path
)
self
.
assertEquals
(
radios_urls_evolution
[
0
].
data_slice
(),
(
'0'
,)
)
self
.
assertEquals
(
radios_urls_evolution
[
0
].
data_path
(),
'
/
entry
_
0000/
measurement/pcoedge64
/data'
)
self
.
assertEquals
(
len
(
radios_urls_evolution
),
1
503
)
self
.
assertEquals
(
radios_urls_evolution
[
0
].
file_path
(),
self
.
scan
.
master_file
)
self
.
assertEquals
(
radios_urls_evolution
[
0
].
data_slice
(),
22
)
self
.
assertEquals
(
radios_urls_evolution
[
0
].
data_path
(),
'entry0000/
instrument/detector
/data'
)
def
testDarkRefUtils
(
self
):
self
.
assertEqual
(
self
.
scan
.
tomo_n
,
100
)
self
.
assertEqual
(
self
.
scan
.
pixel_size
[
1
],
6.5
)
self
.
assertEqual
(
self
.
scan
.
tomo_n
,
1500
)
pixel_size
=
self
.
scan
.
pixel_size
self
.
assertTrue
(
pixel_size
is
not
None
)
self
.
assertTrue
(
numpy
.
isclose
(
self
.
scan
.
pixel_size
,
0.05
*
metricsystem
.
MetricSystem
.
MICROMETER
.
value
))
self
.
assertTrue
(
numpy
.
isclose
(
self
.
scan
.
get_pixel_size
(
unit
=
'micrometer'
),
0.05
))
self
.
assertTrue
(
numpy
.
isclose
(
self
.
scan
.
x_pixel_size
,
0.05
*
metricsystem
.
MetricSystem
.
MICROMETER
.
value
))
self
.
assertTrue
(
numpy
.
isclose
(
self
.
scan
.
y_pixel_size
,
0.05
*
metricsystem
.
MetricSystem
.
MICROMETER
.
value
))
def
testNabuUtil
(
self
):
self
.
assertTrue
(
numpy
.
isclose
(
self
.
scan
.
distance
,
-
19.9735
))
self
.
assertTrue
(
numpy
.
isclose
(
self
.
scan
.
get_distance
(
unit
=
'cm'
),
-
1997.35
))
def
suite
():
...
...
tomoscan/scanbase.py
View file @
2bc2c976
...
...
@@ -32,6 +32,7 @@ import os
import
logging
from
typing
import
Union
from
collections
import
OrderedDict
from
.unitsystem.metricsystem
import
MetricSystem
logger
=
logging
.
getLogger
(
__name__
)
...
...
@@ -45,9 +46,9 @@ class TomoScanBase:
:param scan: path to the root folder containing the scan.
:type: Union[str,None]
"""
_
DICT_TYPE_KEY
=
'type'
DICT_TYPE_KEY
=
'type'
_
DICT_PATH_KEY
=
'path'
DICT_PATH_KEY
=
'path'
_SCHEME
=
None
"""scheme to read data url for this type of acquisition"""
...
...
@@ -56,6 +57,11 @@ class TomoScanBase:
self
.
path
=
scan
self
.
_type
=
type_
def
clear_caches
(
self
):
"""clear caches. Might be call if some data changed after
first read of data or metadata"""
raise
NotImplementedError
(
'Base class'
)
@
property
def
path
(
self
)
->
Union
[
None
,
str
]:
"""
...
...
@@ -111,7 +117,7 @@ class TomoScanBase:
self
.
_flats
=
flats
@
property
def
darks
(
self
)
->
Union
[
None
,
dict
]:
def
darks
(
self
)
->
Union
[
None
,
dict
]:
"""list of darks files"""
return
self
.
_darks
...
...
@@ -134,6 +140,7 @@ class TomoScanBase:
@
property
def
tomo_n
(
self
)
->
Union
[
None
,
int
]:
"""number of projection WITHOUT the return projections"""
raise
NotImplementedError
(
'Base class'
)
@
property
...
...
@@ -141,9 +148,15 @@ class TomoScanBase:
raise
NotImplementedError
(
'Base class'
)
@
property
def
pixel_size
(
self
)
->
Union
[
None
,
in
t
]:
def
pixel_size
(
self
)
->
Union
[
None
,
floa
t
]:
raise
NotImplementedError
(
'Base class'
)
def
get_pixel_size
(
self
,
unit
=
'm'
)
->
Union
[
None
,
float
]:
if
self
.
pixel_size
:
return
self
.
pixel_size
/
MetricSystem
.
from_value
(
unit
).
value
else
:
return
None
@
property
def
dim_1
(
self
)
->
Union
[
None
,
int
]:
raise
NotImplementedError
(
'Base class'
)
...
...
@@ -160,6 +173,33 @@ class TomoScanBase:
def
scan_range
(
self
)
->
Union
[
None
,
int
]:
raise
NotImplementedError
(
'Base class'
)
@
property
def
energy
(
self
)
->
Union
[
None
,
float
]:
"""
:return: incident beam energy in keV
"""
raise
NotImplementedError
(
'Base class'
)
@
property
def
distance
(
self
)
->
Union
[
None
,
float
]:
"""
:return: sample / detector distance in meter
"""
raise
NotImplementedError
(
'Base class'
)
def
get_distance
(
self
,
unit
=
'm'
)
->
Union
[
None
,
float
]:
"""
:param Union[MetricSystem, str] unit: unit requested for the distance
:return: sample / detector distance with the requested unit
"""
if
self
.
distance
:
return
self
.
distance
/
MetricSystem
.
from_value
(
unit
).
value
else
:
return
None
def
update
(
self
)
->
None
:
"""Parse the root folder and files to update informations"""
raise
NotImplementedError
(
"Base class"
)
...
...
@@ -172,8 +212,8 @@ class TomoScanBase:
:rtype: dict
"""
res
=
dict
()
res
[
self
.
_
DICT_TYPE_KEY
]
=
self
.
type
res
[
self
.
_
DICT_PATH_KEY
]
=
self
.
path
res
[
self
.
DICT_TYPE_KEY
]
=
self
.
type
res
[
self
.
DICT_PATH_KEY
]
=
self
.
path
return
res
def
load_from_dict
(
self
,
_dict
:
dict
):
...
...
@@ -203,9 +243,11 @@ class TomoScanBase:
def
get_proj_angle_url
(
self
)
->
dict
:
"""
Return the 'extra' radios of a scan which are used to see if the scan
moved during the acquisition. If no extra radio are found, return the
dictionary of all radios.
return a dictionary of all the projection. key is the angle of the
projection and value is the url.
Keys are int for 'standard' projections and strings for return
projections.
:return dict: angles as keys, radios as value.
"""
...
...
@@ -226,10 +268,10 @@ class TomoScanBase:
:warning: each url should contain only one radio.
:param urls:
ordered lis
t with all the urls. First url should be
:param urls:
dic
t with all the urls. First url should be
the first radio acquire, last url should match the last
radio acquire.
:type:
lis
t
:type:
dic
t
:param n_projection: number of projection for the sample.
:type: int
:param scan_range: acquisition range (usually 180 or 360)
...
...
@@ -241,16 +283,18 @@ class TomoScanBase:
are incoherent
"""
assert
n_projection
is
not
None
ordered_url
=
OrderedDict
(
sorted
(
urls
.
items
(),
key
=
lambda
x
:
x
))
ordered_url
=
OrderedDict
(
sorted
(
urls
.
items
(),
key
=
lambda
x
:
x
))
res
=
{}
# deal with the 'standard' acquisitions
for
proj_i
in
range
(
n_projection
):
angle
=
proj_i
*
scan_range
/
(
n_projection
-
1
)
url
=
list
(
ordered_url
.
values
())[
proj_i
]
if
n_projection
==
1
:
angle
=
0.0
else
:
angle
=
proj_i
*
scan_range
/
(
n_projection
-
1
)
if
proj_i
<
len
(
urls
):
res
[
angle
]
=
urls
[
proj_i
]
# TODO: better to have them from the header ?!!
# haven't found the information in here.
res
[
angle
]
=
url
if
len
(
urls
)
>
n_projection
:
# deal with extra images (used to check if the sampled as moved for
...
...
@@ -261,27 +305,27 @@ class TomoScanBase:
logger
.
warning
(
'incoherent data information to retrieve'
'scan extra images angle'
)
elif
len
(
extraImgs
)
==
4
:
res
[
'270(1)'
]
=
extraImgs
[
0
]
res
[
'180(1)'
]
=
extraImgs
[
1
]
res
[
'90(1)'
]
=
extraImgs
[
2
]
res
[
'0(1)'
]
=
extraImgs
[
3
]
res
[
'270(1)'
]
=
ordered_url
[
extraImgs
[
0
]
]
res
[
'180(1)'
]
=
ordered_url
[
extraImgs
[
1
]
]
res
[
'90(1)'
]
=
ordered_url
[
extraImgs
[
2
]
]
res
[
'0(1)'
]
=
ordered_url
[
extraImgs
[
3
]
]
else
:
res
[
'360(1)'
]
=
extraImgs
[
0
]
res
[
'270(1)'
]
=
extraImgs
[
1
]
res
[
'180(1)'
]
=
extraImgs
[
2
]
res
[
'90(1)'
]
=
extraImgs
[
3
]
res
[
'0(1)'
]
=
extraImgs
[
4
]
res
[
'360(1)'
]
=
ordered_url
[
extraImgs
[
0
]
]
res
[
'270(1)'
]
=
ordered_url
[
extraImgs
[
1
]
]
res
[
'180(1)'
]
=
ordered_url
[
extraImgs
[
2
]
]
res
[
'90(1)'
]
=
ordered_url
[
extraImgs
[
3
]
]
res
[
'0(1)'
]
=
ordered_url
[
extraImgs
[
4
]
]
elif
len
(
extraImgs
)
in
(
2
,
3
):
if
scan_range
>
180
:
logger
.
warning
(
'incoherent data information to retrieve'
'scan extra images angle'
)
elif
len
(
extraImgs
)
is
3
:
res
[
'180(1)'
]
=
extraImgs
[
0
]
res
[
'90(1)'
]
=
extraImgs
[
1
]
res
[
'0(1)'
]
=
extraImgs
[
2
]
res
[
'180(1)'
]
=
ordered_url
[
extraImgs
[
0
]
]
res
[
'90(1)'
]
=
ordered_url
[
extraImgs
[
1
]
]
res
[
'0(1)'
]
=
ordered_url
[
extraImgs
[
2
]
]
else
:
res
[
'90(1)'
]
=
extraImgs
[
0
]
res
[
'0(1)'
]
=
extraImgs
[
1
]
res
[
'90(1)'
]
=
ordered_url
[
extraImgs
[
0
]
]
res
[
'0(1)'
]
=
ordered_url
[
extraImgs
[
1
]
]
else
:
raise
ValueError
(
'incoherent data information to retrieve scan'
'extra images angle'
)
...
...
tomoscan/scanfactory.py
View file @
2bc2c976
...
...
@@ -41,11 +41,10 @@ class ScanFactory:
"""