def pct_mumap( datain, scanner_params, hst=[], t0=0, t1=0, itr=2, petopt='ac', faff='', fpet='', fcomment='', outpath='', store_npy = False, store = False, verbose=True ): ''' GET THE MU-MAP from pCT IMAGE (which is in T1w space) * the mu-map will be registered to PET which will be reconstructed for time frame t0-t1 * it f0 and t1 are not given the whole LM dataset will be reconstructed * the reconstructed PET can be attenuation and scatter corrected or NOT using petopt ''' if not os.path.isfile(faff): from niftypet.nipet.prj import mmrrec # histogram the list data if needed if not hst: from niftypet.nipet.lm import mmrhist hst = mmrhist.mmrhist(datain, scanner_params, t0=t0, t1=t1) # constants, transaxial and axial LUTs are extracted Cnt = scanner_params['Cnt'] txLUT = scanner_params['txLUT'] axLUT = scanner_params['axLUT'] # get hardware mu-map if 'hmumap' in datain and os.path.isfile(datain['hmumap']): muh, _, _ = np.load(datain['hmumap'], allow_pickle=True) if verbose: print 'i> loaded hardware mu-map from file:', datain['hmumap'] elif outpath!='': hmupath = os.path.join( os.path.join(outpath,'mumap-hdw'), 'hmumap.npy') if os.path.isfile( hmupath ): muh, _, _ = np.load(hmupath, allow_pickle=True) datain['hmumap'] = hmupath else: raise IOError('Invalid path to the hardware mu-map') else: print 'e> obtain the hardware mu-map first.' raise IOError('Could not find the hardware mu-map. Have you run the routine for hardware mu-map?') if not 'MRT1W#' in datain and not 'T1nii' in datain and not 'T1bc' in datain: print 'e> no MR T1w images required for co-registration!' raise IOError('Missing MR data') # ---------------------------------- # output dictionary mu_dct = {} if not os.path.isfile(faff): # first recon pet to get the T1 aligned to it if petopt=='qnt': # --------------------------------------------- # OPTION 1 (quantitative recon with all corrections using MR-based mu-map) # get UTE object mu-map (may not be in register with the PET data) mudic = obj_mumap(datain, Cnt) muo = mudic['im'] # reconstruct PET image with UTE mu-map to which co-register T1w recout = mmrrec.osemone( datain, [muh, muo], hst, scanner_params, recmod=3, itr=itr, fwhm=0., fcomment=fcomment+'_qntUTE', outpath=os.path.join(outpath, 'PET', 'positioning'), store_img=True) elif petopt=='nac': # --------------------------------------------- # OPTION 2 (recon without any corrections for scatter and attenuation) # reconstruct PET image with UTE mu-map to which co-register T1w muo = np.zeros(muh.shape, dtype=muh.dtype) recout = mmrrec.osemone( datain, [muh, muo], hst, scanner_params, recmod=1, itr=itr, fwhm=0., fcomment=fcomment+'_NAC', outpath=os.path.join(outpath, 'PET', 'positioning'), store_img=True) elif petopt=='ac': # --------------------------------------------- # OPTION 3 (recon with attenuation correction only but no scatter) # reconstruct PET image with UTE mu-map to which co-register T1w mudic = obj_mumap(datain, Cnt, outpath=outpath) muo = mudic['im'] recout = mmrrec.osemone( datain, [muh, muo], hst, scanner_params, recmod=1, itr=itr, fwhm=0., fcomment=fcomment+'_AC', outpath=os.path.join(outpath, 'PET', 'positioning'), store_img=True) fpet = recout.fpet mu_dct['fpet'] = fpet #------------------------------ # get the affine transformation ft1w = nimpa.pick_t1w(datain) try: regdct = nimpa.coreg_spm( fpet, ft1w, outpath=os.path.join(outpath,'PET', 'positioning') ) except: regdct = nimpa.affine_niftyreg( fpet, ft1w, outpath=os.path.join(outpath,'PET', 'positioning'), #fcomment=fcomment, executable = Cnt['REGPATH'], omp = multiprocessing.cpu_count()/2, rigOnly = True, affDirect = False, maxit=5, speed=True, pi=50, pv=50, smof=0, smor=0, rmsk=True, fmsk=True, rfwhm=15., #millilitres rthrsh=0.05, ffwhm = 15., #millilitres fthrsh=0.05, verbose=verbose ) faff = regdct['faff'] #------------------------------ # pCT file name if outpath=='': pctdir = os.path.dirname(datain['pCT']) else: pctdir = os.path.join(outpath, 'mumap-obj') mmraux.create_dir(pctdir) fpct = os.path.join(pctdir, 'pCT_r_tmp'+fcomment+'.nii.gz') #call the resampling routine to get the pCT in place if os.path.isfile( Cnt['RESPATH'] ): cmd = [Cnt['RESPATH'], '-ref', fpet, '-flo', datain['pCT'], '-trans', faff, '-res', fpct, '-pad', '0'] if not verbose: cmd.append('-voff') call(cmd) else: print 'e> path to resampling executable is incorrect!' sys.exit() # get the NIfTI of the pCT nim = nib.load(fpct) A = nim.get_sform() pct = np.float32( nim.get_data() ) pct = pct[:,::-1,::-1] pct = np.transpose(pct, (2, 1, 0)) # convert the HU units to mu-values mu = hu2mu(pct) # get rid of negatives mu[mu<0] = 0 # return image dictionary with the image itself and other parameters mu_dct['im'] = mu mu_dct['affine'] = A mu_dct['faff'] = faff if store: # now save to numpy array and NIfTI in this folder if outpath=='': pctumapdir = os.path.join( datain['corepath'], 'mumap-obj' ) else: pctumapdir = os.path.join(outpath, 'mumap-obj') mmraux.create_dir(pctumapdir) #> Numpy if store_npy: fnp = os.path.join(pctumapdir, 'mumap-pCT.npy') np.save(fnp, (mu, A, fnp)) # numpy # NIfTI fmu = os.path.join(pctumapdir, 'mumap-pCT' +fcomment+ '.nii.gz') nimpa.array2nii(mu[::-1,::-1,:], A, fmu) mu_dct['fim'] = fmu datain['mumapCT'] = fmu return mu_dct
def osemone(datain, mumaps, hst, scanner_params, recmod=3, itr=4, fwhm=0., mask_radius=29., sctsino=np.array([]), outpath='', store_img=False, frmno='', fcomment='', store_itr=[], emmskS=False, ret_sinos=False, attnsino=None, randsino=None, normcomp=None): #---------- sort out OUTPUT ------------ #-output file name for the reconstructed image, initially assume n/a fout = 'n/a' if store_img or store_itr: if outpath == '': opth = os.path.join(datain['corepath'], 'reconstructed') else: opth = outpath mmraux.create_dir(opth) if ret_sinos: return_ssrb = True return_mask = True else: return_ssrb = False return_mask = False #---------- # Get particular scanner parameters: Constants, transaxial and axial LUTs Cnt = scanner_params['Cnt'] txLUT = scanner_params['txLUT'] axLUT = scanner_params['axLUT'] import time from niftypet import nipet # from niftypet.nipet.sct import mmrsct # from niftypet.nipet.prj import mmrhist if Cnt['VERBOSE']: print 'i> reconstruction in mode', recmod # get object and hardware mu-maps muh, muo = mumaps # get the GPU version of the image dims mus = mmrimg.convert2dev(muo + muh, Cnt) if Cnt['SPN'] == 1: snno = Cnt['NSN1'] elif Cnt['SPN'] == 11: snno = Cnt['NSN11'] # remove gaps from the prompt sino psng = mmraux.remgaps(hst['psino'], txLUT, Cnt) #========================================================================= # GET NORM #------------------------------------------------------------------------- if normcomp == None: ncmp, _ = mmrnorm.get_components(datain, Cnt) else: ncmp = normcomp print 'w> using user-defined normalisation components' nsng = mmrnorm.get_sinog(datain, hst, axLUT, txLUT, Cnt, normcomp=ncmp) #========================================================================= #========================================================================= # ATTENUATION FACTORS FOR COMBINED OBJECT AND BED MU-MAP #------------------------------------------------------------------------- #> combine attenuation and norm together depending on reconstruction mode if recmod == 0: asng = np.ones(psng.shape, dtype=np.float32) else: #> check if the attenuation sino is given as an array if isinstance(attnsino, np.ndarray) \ and attnsino.shape==(Cnt['NSN11'], Cnt['NSANGLES'], Cnt['NSBINS']): asng = mmraux.remgaps(attnsino, txLUT, Cnt) print 'i> using provided attenuation factor sinogram' elif isinstance(attnsino, np.ndarray) \ and attnsino.shape==(Cnt['Naw'], Cnt['NSN11']): asng = attnsino print 'i> using provided attenuation factor sinogram' else: asng = np.zeros(psng.shape, dtype=np.float32) petprj.fprj(asng, mus, txLUT, axLUT, np.array([-1], dtype=np.int32), Cnt, 1) #> combine attenuation and normalisation ansng = asng * nsng #========================================================================= #========================================================================= # Randoms #------------------------------------------------------------------------- if isinstance(randsino, np.ndarray): rsino = randsino rsng = mmraux.remgaps(randsino, txLUT, Cnt) else: rsino, snglmap = nipet.randoms(hst, scanner_params) rsng = mmraux.remgaps(rsino, txLUT, Cnt) #========================================================================= #========================================================================= # SCAT #------------------------------------------------------------------------- if recmod == 2: if sctsino.size > 0: ssng = mmraux.remgaps(sctsino, txLUT, Cnt) elif sctsino.size == 0 and os.path.isfile(datain['em_crr']): emd = nimpa.getnii(datain['em_crr']) ssn = nipet.vsm(datain, mumaps, emd['im'], hst, rsino, scanner_params, prcnt_scl=0.1, emmsk=False) ssng = mmraux.remgaps(ssn, txLUT, Cnt) else: print 'e> no emission image available for scatter estimation! check if it' 's present or the path is correct.' sys.exit() else: ssng = np.zeros(rsng.shape, dtype=rsng.dtype) #========================================================================= if Cnt['VERBOSE']: print '\n>------ OSEM (', itr, ') -------\n' #------------------------------------ Sn = 14 # number of subsets #-get one subset to get number of projection bins in a subset Sprj, s = get_subsets14(0, scanner_params) Nprj = len(Sprj) #-init subset array and sensitivity image for a given subset sinoTIdx = np.zeros((Sn, Nprj + 1), dtype=np.int32) #-init sensitivity images for each subset imgsens = np.zeros((Sn, Cnt['SZ_IMY'], Cnt['SZ_IMX'], Cnt['SZ_IMZ']), dtype=np.float32) for n in range(Sn): sinoTIdx[n, 0] = Nprj #first number of projection for the given subset sinoTIdx[n, 1:], s = get_subsets14(n, scanner_params) # sensitivity image petprj.bprj(imgsens[n, :, :, :], ansng[sinoTIdx[n, 1:], :], txLUT, axLUT, sinoTIdx[n, 1:], Cnt) #------------------------------------- #-mask for reconstructed image. anything outside it is set to zero msk = mmrimg.get_cylinder( Cnt, rad=mask_radius, xo=0, yo=0, unival=1, gpu_dim=True) > 0.9 #-init image img = np.ones((Cnt['SZ_IMY'], Cnt['SZ_IMX'], Cnt['SZ_IMZ']), dtype=np.float32) #-decay correction lmbd = np.log(2) / resources.riLUT[Cnt['ISOTOPE']]['thalf'] if Cnt['DCYCRR'] and 't0' in hst and 'dur' in hst: dcycrr = np.exp(lmbd * hst['t0']) * lmbd * hst['dur'] / ( 1 - np.exp(-lmbd * hst['dur'])) # apply quantitative correction to the image qf = ncmp['qf'] / resources.riLUT[Cnt['ISOTOPE']]['BF'] / float( hst['dur']) qf_loc = ncmp['qf_loc'] elif not Cnt['DCYCRR'] and 't0' in hst and 'dur' in hst: dcycrr = 1. # apply quantitative correction to the image qf = ncmp['qf'] / resources.riLUT[Cnt['ISOTOPE']]['BF'] / float( hst['dur']) qf_loc = ncmp['qf_loc'] else: dcycrr = 1. qf = 1. qf_loc = 1. #-affine matrix for the reconstructed images B = mmrimg.image_affine(datain, Cnt) #-time it stime = time.time() # import pdb; pdb.set_trace() #========================================================================= # OSEM RECONSTRUCTION #------------------------------------------------------------------------- for k in trange(itr, disable=not Cnt['VERBOSE'], desc="OSEM"): petprj.osem(img, msk, psng, rsng, ssng, nsng, asng, imgsens, txLUT, axLUT, sinoTIdx, Cnt) if np.nansum(img) < 0.1: print '---------------------------------------------------------------------' print 'w> it seems there is not enough true data to render reasonable image.' print '---------------------------------------------------------------------' #img[:]=0 itr = k break if recmod >= 3 and (((k < itr - 1) and (itr > 1))): # or (itr==1) sct_time = time.time() sct = nipet.vsm(datain, mumaps, mmrimg.convert2e7(img, Cnt), hst, rsino, scanner_params, emmsk=emmskS, return_ssrb=return_ssrb, return_mask=return_mask) if isinstance(sct, dict): ssn = sct['sino'] else: ssn = sct ssng = mmraux.remgaps(ssn, txLUT, Cnt) if Cnt['VERBOSE']: print 'i> scatter time:', (time.time() - sct_time) # save images during reconstruction if requested if store_itr and k in store_itr: im = mmrimg.convert2e7(img * (dcycrr * qf * qf_loc), Cnt) fout = os.path.join(opth, os.path.basename(datain['lm_bf'])[:8] \ + frmno +'_t'+str(hst['t0'])+'-'+str(hst['t1'])+'sec' \ +'_itr'+str(k)+fcomment+'_inrecon.nii.gz') nimpa.array2nii(im[::-1, ::-1, :], B, fout) if Cnt['VERBOSE']: print 'i> recon time:', (time.time() - stime) #========================================================================= if Cnt['VERBOSE']: print 'i> applying decay correction of', dcycrr print 'i> applying quantification factor', qf, 'to the whole image for the frame duration of :', hst[ 'dur'] img *= dcycrr * qf * qf_loc #additional factor for making it quantitative in absolute terms (derived from measurements) #---- save images ----- #-first convert to standard mMR image size im = mmrimg.convert2e7(img, Cnt) #-description text to NIfTI #-attenuation number: if only bed present then it is 0.5 attnum = (1 * (np.sum(muh) > 0.5) + 1 * (np.sum(muo) > 0.5)) / 2. descrip = 'alg=osem'+ \ ';sub=14'+ \ ';att='+str(attnum*(recmod>0))+ \ ';sct='+str(1*(recmod>1))+ \ ';spn='+str(Cnt['SPN'])+ \ ';itr='+str(itr) +\ ';fwhm='+str(fwhm) +\ ';t0='+str(hst['t0']) +\ ';t1='+str(hst['t1']) +\ ';dur='+str(hst['dur']) +\ ';qf='+str(qf) if fwhm > 0: im = ndi.filters.gaussian_filter(im, fwhm2sig(fwhm, Cnt), mode='mirror') if store_img: fout = os.path.join(opth, os.path.basename(datain['lm_bf'])[:8] \ + frmno +'_t'+str(hst['t0'])+'-'+str(hst['t1'])+'sec' \ +'_itr'+str(itr)+fcomment+'.nii.gz') if Cnt['VERBOSE']: print 'i> saving image to: ', fout nimpa.array2nii(im[::-1, ::-1, :], B, fout, descrip=descrip) # returning: # (0) E7 image [can be smoothed]; # (1) file name of saved E7 image # (2) [optional] scatter sino # (3) [optional] single slice rebinned scatter # (4) [optional] mask for scatter scaling based on attenuation data # (5) [optional] random sino # if ret_sinos and recmod>=3: # recout = namedtuple('recout', 'im, fpet, ssn, sssr, amsk, rsn') # recout.im = im # recout.fpet = fout # recout.ssn = ssn # recout.sssr = sssr # recout.amsk = amsk # recout.rsn = rsino # else: # recout = namedtuple('recout', 'im, fpet') # recout.im = im # recout.fpet = fout if ret_sinos and recmod >= 3 and itr > 1: RecOut = namedtuple('RecOut', 'im, fpet, affine, ssn, sssr, amsk, rsn') recout = RecOut(im, fout, B, ssn, sct['ssrb'], sct['mask'], rsino) else: RecOut = namedtuple('RecOut', 'im, fpet, affine') recout = RecOut(im, fout, B) return recout
def get_hmupos(datain, parts, Cnt, outpath=''): # check if registration executable exists if not os.path.isfile(Cnt['RESPATH']): print 'e> no registration executable found!' sys.exit() #----- get positions from the DICOM list-mode file ----- ihdr, csainfo = mmraux.hdr_lm(datain, Cnt) #table position origin fi = csainfo.find('TablePositionOrigin') tpostr = csainfo[fi:fi+200] tpo = re.sub(r'[^a-zA-Z0-9\-\.]', '', tpostr).split('M') tpozyx = np.array([float(tpo[-1]), float(tpo[-2]), float(tpo[-3])]) / 10 if Cnt['VERBOSE']: print 'i> table position (z,y,x) (cm):', tpozyx #-------------------------------------------------------- #------- get positions from the DICOM mu-map file ------- csamu, dhdr = hdr_mu(datain, Cnt) tmp = re.search('GantryTableHomeOffset(?!_)', csamu) gtostr1 = csamu[ tmp.start():tmp.start()+300 ] gtostr2 = re.sub(r'[^a-zA-Z0-9\-\.]', '', gtostr1) # gantry table offset, through conversion of string to float gtoxyz = re.findall(r'(?<=M)-*[\d]{1,4}\.[\d]{6,9}', gtostr2) gtozyx = np.float32(gtoxyz)[::-1]/10 #-------------------------------------------------------- if Cnt['VERBOSE']: print 'i> gantry table offset (z,y,x) (cm):', gtozyx ## ---- ## old II # csamu, dhdr = nipet.img.mmrimg.hdr_mu(datain, Cnt) # tmp = re.search('GantryTableHomeOffset(?!_)', csamu) # gtostr = csamu[ tmp.start():tmp.start()+300 ] # gto = re.sub(r'[^a-zA-Z0-9\-\.]', '', gtostr).split('M') # # get the first three numbers # zyx = np.zeros(3, dtype=np.float32) # c = 0 # for i in range(len(gto)): # if re.search(r'[\d]{1,3}\.[\d]{6}', gto[i])!=None and c<3: # zyx[c] = np.float32(re.sub(r'[^0-9\-\.]', '', gto[i])) # c+=1 # #gantry table offset # gtozyx = zyx[::-1]/10 ## ---- ## ---- ## old I: only worked for syngo MR B20P # fi = csamu.find('GantryTableHomeOffset') # gtostr =csamu[fi:fi+300] # if dhdr[0x0018, 0x1020].value == 'syngo MR B20P': # gto = re.sub(r'[^a-zA-Z0-9\-\.]', '', gtostr).split('M') # # get the first three numbers # zyx = np.zeros(3, dtype=np.float32) # c = 0 # for i in range(len(gto)): # if re.search(r'[\d]', gto[i])!=None and c<3: # zyx[c] = np.float32(re.sub(r'[^0-9\-\.]', '', gto[i])) # c+=1 # #gantry table offset # gtozyx = zyx[::-1]/10 # if Cnt['VERBOSE']: print 'i> gantry table offset (z,y,x) (cm):', gtozyx # # older scanner version # elif dhdr[0x0018, 0x1020].value == 'syngo MR B18P': # zyx = np.zeros(3, dtype=np.float32) # for k in range(3): # tmp = re.search(r'\{\s*[\-0-9.]*\s*\}', gtostr) # i0 = tmp.start() # i1 = tmp.end() # if gtostr[i0+1:i1-1]!=' ': zyx[k] = np.float32(gtostr[i0+1:i1-1]) # gtostr = gtostr[i1:] # #gantry table offset # gtozyx = zyx[::-1]/10 # if Cnt['VERBOSE']: print 'i> gantry table offset (z,y,x) (cm):', gtozyx ## ----- # create the folder for hardware mu-maps if outpath=='': dirhmu = os.path.join( datain['corepath'], 'mumap-hdw') else: dirhmu = os.path.join( outpath, 'mumap-hdw') mmraux.create_dir(dirhmu) # get the reference nii image fref = os.path.join(dirhmu, 'hmuref.nii.gz') #start horizontal bed position p = re.compile(r'start horizontal bed position.*\d{1,3}\.*\d*') m = p.search(ihdr) fi = ihdr[m.start():m.end()].find('=') hbedpos = 0.1*float(ihdr[m.start()+fi+1:m.end()]) #start vertical bed position p = re.compile(r'start vertical bed position.*\d{1,3}\.*\d*') m = p.search(ihdr) fi = ihdr[m.start():m.end()].find('=') vbedpos = 0.1*float(ihdr[m.start()+fi+1:m.end()]) if Cnt['VERBOSE']: print 'i> creating reference nii image for resampling' B = np.diag(np.array([-10*Cnt['SO_VXX'], 10*Cnt['SO_VXY'], 10*Cnt['SO_VXZ'], 1])) B[0,3] = 10*(.5*Cnt['SO_IMX'])*Cnt['SO_VXX'] B[1,3] = 10*( -.5*Cnt['SO_IMY']+1)*Cnt['SO_VXY'] B[2,3] = 10*((-.5*Cnt['SO_IMZ']+1)*Cnt['SO_VXZ'] + hbedpos ) nimpa.array2nii( np.zeros((Cnt['SO_IMZ'], Cnt['SO_IMY'], Cnt['SO_IMX']), dtype=np.float32), B, fref) #define a dictionary of all positions/offsets of hardware mu-maps hmupos = [None]*5 hmupos[0] = { 'TabPosOrg' : tpozyx, #from DICOM of LM file 'GanTabOff' : gtozyx, #from DICOM of mMR mu-map file 'HBedPos' : hbedpos, #from Interfile of LM file [cm] 'VBedPos' : vbedpos, #from Interfile of LM file [cm] 'niipath' : fref } #-------------------------------------------------------------------------- # iteratively go through the mu-maps and add them as needed for i in parts: fh = os.path.join(Cnt['HMUDIR'], Cnt['HMULIST'][i-1]) # get the interfile header and binary data hdr, im = rd_hmu(fh) #get shape, origin, offset and voxel size s = hmu_shape(hdr) im.shape = s # get the origin, offset and voxel size for the mu-map interfile data org = hmu_origin(hdr) off = hmu_offset(hdr) vs = hmu_voxsize(hdr) # corner voxel position for the interfile image data vpos = (-org*vs + off + gtozyx - tpozyx) #add to the dictionary hmupos[i] = { 'vpos' : vpos, 'shape' : s, #from interfile 'iorg' : org, #from interfile 'ioff' : off, #from interfile 'ivs' : vs, #from interfile 'img' : im, #from interfile 'niipath' : os.path.join(dirhmu, '_'+Cnt['HMULIST'][i-1].split('.')[0]+'.nii.gz') } #save to NIfTI if Cnt['VERBOSE']: print 'i> creating mu-map for:', Cnt['HMULIST'][i-1] A = np.diag(np.append(10*vs[::-1], 1)) A[0,0] *= -1 A[0,3] = 10*(-vpos[2]) A[1,3] = -10*((s[1]-1)*vs[1] + vpos[1]) A[2,3] = -10*((s[0]-1)*vs[0] - vpos[0]) nimpa.array2nii(im[::-1,::-1,:], A, hmupos[i]['niipath']) # resample using nify.reg fout = os.path.join( os.path.dirname (hmupos[0]['niipath']), 'r'+os.path.basename(hmupos[i]['niipath']).split('.')[0]+'.nii.gz' ) cmd = [ Cnt['RESPATH'], '-ref', hmupos[0]['niipath'], '-flo', hmupos[i]['niipath'], '-res', fout, '-pad', '0'] if not Cnt['VERBOSE']: cmd.append('-voff') call(cmd) return hmupos