# -*- coding = utf-8 -*- """ SimNICT Dataset Generator Generates Non-Ideal measurement CT (NICT) simulations from preprocessed ICT data This script creates three types of NICT simulations: 1. Sparse-View CT (SVCT): Limited projection views (15-360 views) 2. Limited-Angle CT (LACT): Restricted angular range (75°-270°) 3. Low-Dose CT (LDCT): Reduced photon dose (5%-75% of normal dose) File Structure: - Input: dataset_name/volume_xxx.nii.gz (flat structure from Internet Archive) - Output: output_path/dataset_name/volume_xxx.nii.gz Usage: python simnict_generator.py Dependencies: numpy, torch, nibabel, odl, astra-toolbox, opencv-python, pillow """ from __future__ import absolute_import, print_function import numpy as np import time import os import nibabel as nib import odl import random import astra # Dataset configuration DATASETS = ['AMOS', 'COVID_19_NY_SBU', 'CT Images in COVID-19', 'CT_COLONOGRAPHY', 'LNDb', 'LUNA', 'MELA', 'STOIC'] # Path configuration INPUT_PATH = 'M:/' # Original ICT data path (downloaded from Internet Archive) OUTPUT_SVCT = 'K:/SpV/' # Sparse-view CT output path OUTPUT_LACT = 'K:/LmV/' # Limited-angle CT output path OUTPUT_LDCT = 'K:/LD/' # Low-dose CT output path # Simulation parameter ranges SVCT_VIEW_RANGE = (15, 360) # Number of projection views LACT_ANGLE_RANGE = (75, 270) # Angular range in degrees LDCT_DOSE_RANGE = (5, 75) # Dose percentage def process_dataset(input_path, dataset_name): """ Process a complete dataset to generate NICT simulations Args: input_path (str): Path to input ICT data dataset_name (str): Name of the dataset to process """ ict_path = os.path.join(input_path, dataset_name) if not os.path.exists(ict_path): print(f"Warning: Path {ict_path} does not exist, skipping {dataset_name}") return # Create output directories for output_path in [OUTPUT_SVCT, OUTPUT_LACT, OUTPUT_LDCT]: os.makedirs(os.path.join(output_path, dataset_name), exist_ok=True) files = os.listdir(ict_path) num_files = len(files) print(f"Processing {dataset_name}: {num_files} files") for i, filename in enumerate(files): print(f"Processing {dataset_name} - File {i+1}/{num_files}: {filename}") # Load ICT volume ict_file_path = os.path.join(ict_path, filename) image_obj = nib.load(ict_file_path) ict_volume = image_obj.get_fdata() + 1024 # Convert to [0, 4096] range ict_volume[ict_volume < 0] = 0 L, W, S = ict_volume.shape # Initialize output volumes svct_volume = np.zeros((L, W, S), dtype=np.int16) lact_volume = np.zeros((L, W, S), dtype=np.int16) ldct_volume = np.zeros((L, W, S), dtype=np.int16) # Process each slice for slice_idx in range(S): ict_slice = ict_volume[:, :, slice_idx] # Generate SVCT with random view number svct_views = random.randint(*SVCT_VIEW_RANGE) svct_slice = create_sparse_view_ct(ict_slice, L, W, svct_views) svct_volume[:, :, slice_idx] = svct_slice # Generate LACT with random angular range lact_angle = random.randint(*LACT_ANGLE_RANGE) lact_slice = create_limited_angle_ct(ict_slice, L, W, lact_angle) lact_volume[:, :, slice_idx] = lact_slice # Generate LDCT with random dose level ldct_dose = random.randint(*LDCT_DOSE_RANGE) ldct_slice = create_low_dose_ct(ict_slice - 1024, L, W, ldct_dose) # Convert back to [-1024, 3072] ldct_volume[:, :, slice_idx] = ldct_slice # Save NICT volumes save_nict_volume(svct_volume, OUTPUT_SVCT, dataset_name, filename, volume_type="SVCT") save_nict_volume(lact_volume, OUTPUT_LACT, dataset_name, filename, volume_type="LACT") save_nict_volume(ldct_volume, OUTPUT_LDCT, dataset_name, filename, volume_type="LDCT") def save_nict_volume(volume, output_path, dataset_name, filename, volume_type): """Save NICT volume to NIfTI format""" if volume_type in ["SVCT", "LACT"]: # Convert from [0, 4096] to [-1024, 3072] range volume_output = volume - 1024 volume_output[volume_output < -1024] = -1024 else: # LDCT # Already in [-1024, 3072] range volume_output = volume volume_output[volume_output < -1024] = -1024 nifti_image = nib.Nifti1Image(volume_output, np.eye(4)) output_file = os.path.join(output_path, dataset_name, filename) nib.save(nifti_image, output_file) def create_sparse_view_ct(ict_slice, height, width, num_views): """ Generate Sparse-View CT using ODL Args: ict_slice: Input ICT slice [0, 4096] height, width: Image dimensions num_views: Number of projection views Returns: Reconstructed sparse-view CT slice """ # Create reconstruction space reco_space = odl.uniform_discr( min_pt=[-height/4, -width/4], max_pt=[height/4, width/4], shape=[height, width], dtype='float32' ) # Define geometry with limited views angle_partition = odl.uniform_partition(0, 2 * np.pi, num_views) detector_partition = odl.uniform_partition(-360, 360, 1024) geometry = odl.tomo.FanBeamGeometry( angle_partition, detector_partition, src_radius=1270, det_radius=870 ) # Create ray transform and reconstruct ray_trafo = odl.tomo.RayTransform(reco_space, geometry) projection = ray_trafo(ict_slice.astype('float32')) fbp = odl.tomo.fbp_op(ray_trafo) reconstruction = fbp(projection) return reconstruction def create_limited_angle_ct(ict_slice, height, width, angle_range): """ Generate Limited-Angle CT using ODL Args: ict_slice: Input ICT slice [0, 4096] height, width: Image dimensions angle_range: Angular range in degrees Returns: Reconstructed limited-angle CT slice """ # Create reconstruction space reco_space = odl.uniform_discr( min_pt=[-height/4, -width/4], max_pt=[height/4, width/4], shape=[height, width], dtype='float32' ) # Define geometry with limited angular range angle_fraction = angle_range / 360 num_angles = int(720 * angle_fraction) angle_partition = odl.uniform_partition(0, 2 * np.pi * angle_fraction, num_angles) detector_partition = odl.uniform_partition(-360, 360, 1024) geometry = odl.tomo.FanBeamGeometry( angle_partition, detector_partition, src_radius=1270, det_radius=870 ) # Create ray transform and reconstruct ray_trafo = odl.tomo.RayTransform(reco_space, geometry) projection = ray_trafo(ict_slice.astype('float32')) fbp = odl.tomo.fbp_op(ray_trafo) reconstruction = fbp(projection) return reconstruction def create_low_dose_ct(ict_slice, height, width, dose_percentage): """ Generate Low-Dose CT using ASTRA with Poisson noise simulation Args: ict_slice: Input ICT slice [-1024, 3072] height, width: Image dimensions dose_percentage: Dose level as percentage of normal dose Returns: Reconstructed low-dose CT slice """ dose_fraction = dose_percentage / 100.0 # Convert to attenuation coefficients u = 0.0192 # Linear attenuation coefficient attenuation_map = ict_slice * u / 1000.0 + u # ASTRA geometry setup vol_geom = astra.create_vol_geom([height, width]) angles = np.linspace(np.pi, -np.pi, 720) proj_geom = astra.create_proj_geom( 'fanflat', 1.685839319229126, 1024, angles, 600.4500331878662, 485.1499423980713 ) # Create projector and forward project proj_id = astra.create_projector('cuda', proj_geom, vol_geom) operator = astra.OpTomo(proj_id) # Forward projection sinogram = operator * np.mat(attenuation_map) / 2 # Add Poisson noise based on dose level noise = np.random.normal(0, 1, 720 * 1024) noise_scaling = np.sqrt((1 - dose_fraction) / dose_fraction * (np.exp(sinogram) / 1e6)) noisy_sinogram = sinogram + noise * noise_scaling # Reconstruct with FBP noisy_sinogram_2d = np.reshape(noisy_sinogram, [720, -1]) reconstruction = operator.reconstruct('FBP_CUDA', noisy_sinogram_2d) # Convert back to HU values reconstruction = reconstruction.reshape((height, width)) return (reconstruction * 2 - u) / u * 1000 def main(): """Main processing function""" print('SimNICT Dataset Generator Started') start_time = time.time() # Process each dataset for dataset_name in DATASETS: print(f"\n{'='*50}") print(f"Processing Dataset: {dataset_name}") print(f"{'='*50}") try: process_dataset(INPUT_PATH, dataset_name) duration = time.time() - start_time print(f"Completed {dataset_name} in {duration:.1f} seconds") except Exception as e: print(f"Error processing {dataset_name}: {str(e)}") continue total_duration = time.time() - start_time print(f"\nSimNICT Generation Complete - Total time: {total_duration/3600:.2f} hours") if __name__ == "__main__": main()