def compareTwoVolumes(particle,reference,referenceWeighting,wedgeInfo,rotations,shift, scoreObject=0,mask=None,preprocessing=[],binning=1,verbose=False): """ compare: Compares two volumes @param particle: A particle @param reference: A reference @param referenceWeighting: Fourier weighting of the reference (sum of wedges for instance) @param wedgeInfo: What does the wedge look alike? @type wedgeInfo: L{pytom.basic.structures.WedgeInfo} @param rotations: All rotations to be scanned @type rotations: Must be a L{pytom.angles.angleList.OneAngleList} @param scoreObject: @type scoreObject: L{pytom.score.score.Score} @param mask: @param preprocessing: Class storing preprocessing of particle and reference such as bandpass @type preprocessing: L{pytom.alignment.preprocessing.Preprocessing} @param binning: Is binning applied (Disabled when == 1) @return: Returns the correlation coefficient of two volumes. @author: Thomas Hrabe """ from pytom_volume import vol,transformSpline from pytom.basic.filter import filter,rotateWeighting from pytom.angles.angleList import OneAngleList assert particle.__class__ == vol, "particle not of type vol" assert reference.__class__ == vol, "reference not of type vol" assert referenceWeighting.__class__ == vol or referenceWeighting.__class__ == str assert rotations.__class__ == OneAngleList or len(rotations) == 1 if scoreObject == 0: from pytom.score.score import nxcfScore scoreObject = nxcfScore() particleCopy = vol(particle.sizeX(),particle.sizeY(),particle.sizeZ()) #process particle particle = wedgeInfo.apply(particle) particle = preprocessing.apply(particle,True) particleCopy.copyVolume(particle) #manipulate reference currentRotation = rotations.nextRotation() sizeX = particle.sizeX() sizeY = particle.sizeY() sizeZ = particle.sizeZ() centerX = sizeX/2.0 centerY = sizeY/2.0 centerZ = sizeZ/2.0 simulatedVol= vol(sizeX,sizeY,sizeZ) transformSpline(reference,simulatedVol,currentRotation[0],currentRotation[1],currentRotation[2], centerX,centerY,centerZ,0,0,0,shift[0]/binning,shift[1]/binning,shift[2]/binning) simulatedVol = wedgeInfo.apply(simulatedVol) m = mask.getVolume(currentRotation) simulatedVol = simulatedVol * m simulatedVol = preprocessing.apply(simulatedVol,True) if not referenceWeighting.__class__ == str: from pytom_freqweight import weight weightingRotated = rotateWeighting(weighting=referenceWeighting, z1=currentRotation[0], z2=currentRotation[1], x=currentRotation[2], isReducedComplex=None, returnReducedComplex=True, binarize=False) if verbose: particleCopy.info('pc') weightingRotated.info('wr') r = list(filter(particleCopy,weight(weightingRotated))) particleCopy = r[0] scoringResult = scoreObject.scoringCoefficient(particleCopy,simulatedVol,m) if scoringResult.__class__ == list: scoringResult = scoringResult[0] assert scoringResult == scoringResult, "scoringResult possibly NaN" return scoringResult
def average2(particleList, weighting=False, norm=False, determine_resolution=False, mask=None, binning=1, verbose=False): """ 2nd version of average function. Will not write the averages to the disk. Also support internal \ resolution determination. """ from pytom_volume import read, vol, complexDiv, complexRealMult from pytom_volume import transformSpline as transform from pytom.basic.fourier import fft, ifft, convolute from pytom.basic.normalise import mean0std1 from pytom.tools.ProgressBar import FixedProgBar from pytom.basic.filter import lowpassFilter, rotateWeighting from math import exp if len(particleList) == 0: raise RuntimeError('The particlelist provided is empty. Aborting!') if verbose: progressBar = FixedProgBar(0,len(particleList),'Particles averaged ') progressBar.update(0) numberAlignedParticles = 0 even = None odd = None wedgeSum_even = None wedgeSum_odd = None newParticle = None is_odd = True for particleObject in particleList: particle = read(particleObject.getFilename(), 0,0,0,0,0,0,0,0,0, binning,binning,binning) if norm: mean0std1(particle) wedgeInfo = particleObject.getWedge() # apply its wedge to itself particle = wedgeInfo.apply(particle) if odd is None: # initialization sizeX = particle.sizeX() sizeY = particle.sizeY() sizeZ = particle.sizeZ() newParticle = vol(sizeX,sizeY,sizeZ) centerX = sizeX/2 centerY = sizeY/2 centerZ = sizeZ/2 odd = vol(sizeX,sizeY,sizeZ) odd.setAll(0.0) even = vol(sizeX,sizeY,sizeZ) even.setAll(0.0) wedgeSum_odd = wedgeInfo.returnWedgeVolume(sizeX,sizeY,sizeZ) wedgeSum_odd.setAll(0) wedgeSum_even = wedgeInfo.returnWedgeVolume(sizeX,sizeY,sizeZ) wedgeSum_even.setAll(0) # create spectral wedge weighting rotation = particleObject.getRotation() rotinvert = rotation.invert() if analytWedge: # > original buggy version wedge = wedgeInfo.returnWedgeVolume(sizeX,sizeY,sizeZ,False, rotinvert) # < original buggy version else: # > FF: interpol bugfix wedge = rotateWeighting( weighting=wedgeInfo.returnWedgeVolume(sizeX,sizeY,sizeZ,False), z1=rotinvert[0], z2=rotinvert[1], x=rotinvert[2], mask=None, isReducedComplex=True, returnReducedComplex=True) # < FF # > TH bugfix #wedgeVolume = wedgeInfo.returnWedgeVolume(wedgeSizeX=sizeX, wedgeSizeY=sizeY, wedgeSizeZ=sizeZ, # humanUnderstandable=True, rotation=rotinvert) #wedge = rotate(volume=wedgeVolume, rotation=rotinvert, imethod='linear') # < TH if is_odd: wedgeSum_odd = wedgeSum_odd + wedge else: wedgeSum_even = wedgeSum_even + wedge # shift and rotate particle shiftV = particleObject.getShift() newParticle.setAll(0) transform(particle,newParticle,-rotation[1],-rotation[0],-rotation[2], centerX,centerY,centerZ,-shiftV[0]/binning, -shiftV[1]/binning,-shiftV[2]/binning,0,0,0) if is_odd: if weighting: weight = 1. - particleObject.getScore().getValue() #weight = weight**2 weight = exp(-1.*weight) odd = odd + newParticle * weight else: odd = odd + newParticle else: if weighting: weight = 1. - particleObject.getScore().getValue() #weight = weight**2 weight = exp(-1.*weight) even = even + newParticle * weight else: even = even + newParticle is_odd = not is_odd if verbose: numberAlignedParticles = numberAlignedParticles + 1 progressBar.update(numberAlignedParticles) # determine resolution if needed fsc = None if determine_resolution: # apply spectral weighting to sum f_even = fft(even) w_even = complexDiv(f_even, wedgeSum_even) w_even = ifft(w_even) w_even.shiftscale(0.0,1/float(sizeX*sizeY*sizeZ)) f_odd = fft(odd) w_odd = complexDiv(f_odd, wedgeSum_odd) w_odd = ifft(w_odd) w_odd.shiftscale(0.0,1/float(sizeX*sizeY*sizeZ)) from pytom.basic.correlation import FSC fsc = FSC(w_even, w_odd, sizeX/2, mask, verbose=False) # add together result = even+odd wedgeSum = wedgeSum_even+wedgeSum_odd invert_WedgeSum( invol=wedgeSum, r_max=sizeX/2-2., lowlimit=.05*len(particleList), lowval=.05*len(particleList)) #wedgeSum.write(averageName[:len(averageName)-3] + '-WedgeSumInverted.em') result = convolute(v=result, k=wedgeSum, kernel_in_fourier=True) # do a low pass filter #result = lowpassFilter(result, sizeX/2-2, (sizeX/2-1)/10.)[0] return (result, fsc)
def bestAlignment(particle, reference, referenceWeighting, wedgeInfo, rotations, scoreObject=0, mask=None, preprocessing=None, progressBar=False, binning=1, bestPeak=None, verbose=False): """ bestAlignment: Determines best alignment of particle relative to the reference @param particle: A particle @type particle: L{pytom_volume.vol} @param reference: A reference @type reference: L{pytom_volume.vol} @param referenceWeighting: Fourier weighting of the reference (sum of wedges for instance) @type referenceWeighting: L{pytom.basic.structures.vol} @param wedgeInfo: What does the wedge look alike? @type wedgeInfo: L{pytom.basic.structures.Wedge} @param rotations: All rotations to be scanned @type rotations: L{pytom.angles.AngleObject} @param scoreObject: @type scoreObject: L{pytom.score.score.Score} @param mask: real-space mask for correlation function @type mask: L{pytom.basic.structures.Particle} @param preprocessing: Class storing preprocessing of particle and reference such as bandpass @type preprocessing: L{pytom.alignment.preprocessing.Preprocessing} @param progressBar: Display progress bar of alignment. False by default. @param binning: binning factor - e.g. binning=2 reduces size by FACTOR of 2 @type binning: int or float @param bestPeak: Initialise best peak with old values. @param verbose: Print out infos. Writes CC volume to disk!!! Default is False @return: Returns the best rotation for particle and the corresponding scoring result. @author: Thomas Hrabe """ from pytom.basic.correlation import subPixelPeak, subPixelPeakParabolic from pytom.alignment.structures import Peak from pytom_volume import peak, vol, vol_comp from pytom.basic.filter import filter,rotateWeighting from pytom.basic.structures import Rotation, Shift, Particle, Mask from pytom.angles.angle import AngleObject from pytom.alignment.preprocessing import Preprocessing from pytom.basic.transformations import resize, resizeFourier binningType = 'Fourier' # or 'Fourier' assert isinstance(rotations, AngleObject), "bestAlignment: rotations must be " \ "AngleObject!" currentRotation = rotations.nextRotation() if currentRotation == [None,None,None]: raise Exception('bestAlignment: No rotations are sampled! Something is wrong with input rotations') assert particle.__class__ == vol, "particle not of type vol" assert reference.__class__ == vol, "reference not of type vol" assert (referenceWeighting.__class__ == vol or referenceWeighting.__class__ == str), \ "referenceWeighting not volume or str" if mask: assert mask.__class__ == Mask, "Mask not of type Mask" m = mask.getVolume() if scoreObject == 0 or not scoreObject: from pytom.score.score import xcfScore scoreObject = xcfScore() # fix binning if binning == 0: binning = 1 if binning != 1: particleUnbinned = vol(particle.sizeX(), particle.sizeY(), particle.sizeZ()) particleUnbinned.copyVolume(particle) particle = resize(volume=particle, factor=1./binning, interpolation=binningType) if type(particle) == tuple: particle = particle[0] referenceUnbinned = vol(reference.sizeX(), reference.sizeY(), reference.sizeZ()) referenceUnbinned.copyVolume(reference) reference = resize(volume=reference, factor=1./binning, interpolation=binningType) if type(reference) == tuple: reference = reference[0] if mask: m = resize(volume=m, factor=1./binning, interpolation='Spline') if not referenceWeighting.__class__ == str: referenceWeightingUnbinned = vol_comp(referenceWeighting.sizeX(), referenceWeighting.sizeY(), referenceWeighting.sizeZ()) referenceWeightingUnbinned.copyVolume(referenceWeighting) if binning != 1: referenceWeighting = resizeFourier(fvol=referenceWeighting, factor=1./binning) centerX, centerY, centerZ = int(particle.sizeX()/2), int(particle.sizeY()/2), int(particle.sizeZ()/2) # create buffer volume for transformed particle particleCopy = vol(particle.sizeX(),particle.sizeY(),particle.sizeZ()) particle = wedgeInfo.apply(particle) #apply wedge to itself if preprocessing is None: preprocessing = Preprocessing() preprocessing.setTaper( taper=particle.sizeX()/10.) particle = preprocessing.apply(volume=particle, bypassFlag=True) # filter particle to some resolution particleCopy.copyVolume(particle) # compute standard veviation volume really only if needed if mask and (scoreObject._type=='FLCFScore'): from pytom_volume import sum from pytom.basic.correlation import meanUnderMask, stdUnderMask p = sum(m) meanV = meanUnderMask(particle, m, p) stdV = stdUnderMask(particle, m, p, meanV) else: meanV = None stdV = None while currentRotation != [None,None,None]: if mask: m = mask.getVolume(currentRotation) if binning != 1: m = resize(volume=m, factor=1./binning, interpolation='Spline') #update stdV if mask is not a sphere # compute standard deviation volume really only if needed if not mask.isSphere() and (scoreObject._type=='FLCFScore'): meanV = meanUnderMask(particle, m, p) stdV = stdUnderMask(particle, m, p, meanV) else: m = None simulatedVol = _rotateWedgeReference(reference, currentRotation, wedgeInfo, m, [centerX, centerY, centerZ]) simulatedVol = preprocessing.apply(volume=simulatedVol, bypassFlag=True) #weight particle if not referenceWeighting.__class__ == str: from pytom_freqweight import weight weightingRotated = rotateWeighting(weighting=referenceWeighting, z1=currentRotation[0], z2=currentRotation[1], x=currentRotation[2], isReducedComplex=True, returnReducedComplex=True, binarize=False) particleCopy.copyVolume(particle) r = list(filter(particleCopy, weight(weightingRotated))) particleCopy = r[0] scoringResult = scoreObject.score(particleCopy, simulatedVol, m, stdV) pk = peak(scoringResult) #with subPixelPeak [peakValue,peakPosition] = subPixelPeak(scoreVolume=scoringResult, coordinates=pk, interpolation='Quadratic', verbose=False) #[peakValue,peakPosition] = subPixelPeakParabolic(scoreVolume=scoringResult, coordinates=pk, verbose=False) # determine shift relative to center shiftX = (peakPosition[0] - centerX) * binning shiftY = (peakPosition[1] - centerY) * binning shiftZ = (peakPosition[2] - centerZ) * binning #NANs would fail this test. assert peakValue == peakValue, "peakValue seems to be NaN" newPeak = Peak(peakValue, Rotation(currentRotation), Shift(shiftX, shiftY, shiftZ)) if verbose: print('Rotation: z1=%3.1f, z2=%3.1f, x=%3.1f; Dx=%2.2f, Dy=%2.2f, Dz=%2.2f, CC=%2.3f' % \ (currentRotation[0], currentRotation[1], currentRotation[2], shiftX, shiftY, shiftZ, peakValue)) if bestPeak is None: bestPeak = newPeak if verbose: scoringResult.write('BestScore.em') if bestPeak < newPeak: bestPeak = newPeak if verbose: scoringResult.write('BestScore.em') currentRotation = rotations.nextRotation() # repeat ccf for binned sampling to get translation accurately if binning != 1: m = mask.getVolume(bestPeak.getRotation()) centerX, centerY, centerZ = int(particleUnbinned.sizeX()/2), int(particleUnbinned.sizeY()/2), \ int(particleUnbinned.sizeZ()/2) simulatedVol = _rotateWedgeReference(referenceUnbinned, bestPeak.getRotation(), wedgeInfo, m, [centerX, centerY, centerZ]) simulatedVol = preprocessing.apply(volume=simulatedVol, bypassFlag=True) if mask and scoreObject._type=='FLCFScore': p = sum(m) meanV = meanUnderMask(volume=particleUnbinned, mask=m, p=p) stdV = stdUnderMask(volume=particleUnbinned, mask=m, p=p, meanV=meanV) scoreObject._peakPrior.reset_weight() scoringResult = scoreObject.score(particle=particleUnbinned, reference=simulatedVol, mask=m, stdV=stdV) pk = peak(scoringResult) [peakValue,peakPosition] = subPixelPeak(scoreVolume=scoringResult, coordinates=pk, interpolation='Quadratic', verbose=False) shiftX = (peakPosition[0] - centerX) shiftY = (peakPosition[1] - centerY) shiftZ = (peakPosition[2] - centerZ) bestPeak = Peak(peakValue, bestPeak.getRotation(), Shift(shiftX, shiftY, shiftZ)) if verbose: print('BestAlignment: z1=%3.1f, z2=%3.1f, x=%3.1f; Dx=%2.2f, Dy=%2.2f, Dz=%2.2f, CC=%2.3f' % \ (bestPeak.getRotation()[0], bestPeak.getRotation()[1], bestPeak.getRotation()[2], bestPeak.getShift()[0], bestPeak.getShift()[1], bestPeak.getShift()[2], bestPeak.getScoreValue())) rotations.reset() scoreObject._peakPrior.reset_weight() return bestPeak
def average( particleList, averageName, showProgressBar=False, verbose=False, createInfoVolumes=False, weighting=False, norm=False): """ average : Creates new average from a particleList @param particleList: The particles @param averageName: Filename of new average @param verbose: Prints particle information. Disabled by default. @param createInfoVolumes: Create info data (wedge sum, inverted density) too? False by default. @param weighting: apply weighting to each average according to its correlation score @param norm: apply normalization for each particle @return: A new Reference object @rtype: L{pytom.basic.structures.Reference} @author: Thomas Hrabe @change: limit for wedgeSum set to 1% or particles to avoid division by small numbers - FF """ from pytom_volume import read,vol,reducedToFull,limit, complexRealMult from pytom.basic.filter import lowpassFilter, rotateWeighting from pytom_volume import transformSpline as transform from pytom.basic.fourier import convolute from pytom.basic.structures import Reference from pytom.basic.normalise import mean0std1 from pytom.tools.ProgressBar import FixedProgBar from math import exp import os if len(particleList) == 0: raise RuntimeError('The particle list is empty. Aborting!') if showProgressBar: progressBar = FixedProgBar(0,len(particleList),'Particles averaged ') progressBar.update(0) numberAlignedParticles = 0 result = [] wedgeSum = [] newParticle = None # pre-check that scores != 0 if weighting: wsum = 0. for particleObject in particleList: wsum += particleObject.getScore().getValue() if wsum < 0.00001: weighting = False print("Warning: all scores have been zero - weighting not applied") for particleObject in particleList: if verbose: print(particleObject) if not os.path.exists(particleObject.getFilename()): continue particle = read(particleObject.getFilename()) if norm: # normalize the particle mean0std1(particle) # happen inplace wedgeInfo = particleObject.getWedge() # apply its wedge to itself particle = wedgeInfo.apply(particle) if result == []: sizeX = particle.sizeX() sizeY = particle.sizeY() sizeZ = particle.sizeZ() newParticle = vol(sizeX,sizeY,sizeZ) centerX = sizeX/2 centerY = sizeY/2 centerZ = sizeZ/2 result = vol(sizeX,sizeY,sizeZ) result.setAll(0.0) if analytWedge: wedgeSum = wedgeInfo.returnWedgeVolume(wedgeSizeX=sizeX, wedgeSizeY=sizeY, wedgeSizeZ=sizeZ) else: # > FF bugfix wedgeSum = wedgeInfo.returnWedgeVolume(sizeX,sizeY,sizeZ) # < FF # > TH bugfix #wedgeSum = vol(sizeX,sizeY,sizeZ) # < TH #wedgeSum.setAll(0) assert wedgeSum.sizeX() == sizeX and wedgeSum.sizeY() == sizeY and wedgeSum.sizeZ() == sizeZ/2+1, \ "wedge initialization result in wrong dims :(" wedgeSum.setAll(0) ### create spectral wedge weighting rotation = particleObject.getRotation() rotinvert = rotation.invert() if analytWedge: # > analytical buggy version wedge = wedgeInfo.returnWedgeVolume(sizeX,sizeY,sizeZ,False, rotinvert) else: # > FF: interpol bugfix wedge = rotateWeighting( weighting=wedgeInfo.returnWedgeVolume(sizeX,sizeY,sizeZ,False), z1=rotinvert[0], z2=rotinvert[1], x=rotinvert[2], mask=None, isReducedComplex=True, returnReducedComplex=True) # < FF # > TH bugfix #wedgeVolume = wedgeInfo.returnWedgeVolume(wedgeSizeX=sizeX, wedgeSizeY=sizeY, wedgeSizeZ=sizeZ, # humanUnderstandable=True, rotation=rotinvert) #wedge = rotate(volume=wedgeVolume, rotation=rotinvert, imethod='linear') # < TH ### shift and rotate particle shiftV = particleObject.getShift() newParticle.setAll(0) transform(particle,newParticle,-rotation[1],-rotation[0],-rotation[2], centerX,centerY,centerZ,-shiftV[0],-shiftV[1],-shiftV[2],0,0,0) if weighting: weight = 1.-particleObject.getScore().getValue() #weight = weight**2 weight = exp(-1.*weight) result = result + newParticle * weight wedgeSum = wedgeSum + wedge * weight else: result = result + newParticle wedgeSum = wedgeSum + wedge if showProgressBar: numberAlignedParticles = numberAlignedParticles + 1 progressBar.update(numberAlignedParticles) ###apply spectral weighting to sum result = lowpassFilter(result, sizeX/2-1, 0.)[0] #if createInfoVolumes: result.write(averageName[:len(averageName)-3]+'-PreWedge.em') wedgeSum.write(averageName[:len(averageName)-3] + '-WedgeSumUnscaled.em') invert_WedgeSum( invol=wedgeSum, r_max=sizeX/2-2., lowlimit=.05*len(particleList), lowval=.05*len(particleList)) if createInfoVolumes: wedgeSum.write(averageName[:len(averageName)-3] + '-WedgeSumInverted.em') result = convolute(v=result, k=wedgeSum, kernel_in_fourier=True) # do a low pass filter #result = lowpassFilter(result, sizeX/2-2, (sizeX/2-1)/10.)[0] result.write(averageName) if createInfoVolumes: resultINV = result * -1 #write sign inverted result to disk (good for chimera viewing ... ) resultINV.write(averageName[:len(averageName)-3]+'-INV.em') newReference = Reference(averageName,particleList) return newReference