Borehole NMR processing
You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

MRProc.py 36KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830
  1. from __future__ import division # in case this is called from python2
  2. from __future__ import print_function
  3. import numpy as np
  4. import matplotlib.pyplot as plt
  5. import sys, os
  6. import notch
  7. from scipy import signal
  8. from scipy.interpolate import splprep, splev, spline # for smoothing spline FITPACK
  9. from scipy.interpolate import UnivariateSpline
  10. from decay import *
  11. import scipy.stats as stats
  12. from pwctimeWhite import *
  13. # for signal slot communication with GUI
  14. #from PyQt4.QtCore import QObject, SIGNAL
  15. from PySide.QtCore import QObject, SIGNAL
  16. from pandas.io.parsers import read_csv #<<- faster than Numpy but not working right now
  17. from timeit import timeit
  18. class MRProc(QObject):
  19. """ Class to process read and invert ORS borehole data
  20. """
  21. def __init__(self):
  22. QObject.__init__(self)
  23. self.burst = False # For Schlumberger or VC burst data
  24. #def __init__(self, filename):
  25. # # what do we need to know that isn't in new files?
  26. # fileName, fileExtension = os.path.splitext(filename)
  27. # if fileExtension == ".ors":
  28. # self.loadORSFile(filename)
  29. #self.DistRes = { }
  30. #self.DistRes = { 'env':np.array(0), 'mod':np.array(0), 'Time':np.array(0) } #[a0,b0,rt20] }
  31. def loadORSFile_OLD(self, filename):
  32. """ Loads ORS files with minimal header info
  33. """
  34. f = open(filename)
  35. header = f.readline().split() #np.loadtxt(filename, comment="#", )
  36. ih = 1
  37. while (header[0] == "#"):
  38. ih += 1
  39. header = f.readline().split() #np.loadtxt(filename, comment="#", )
  40. self.NR = eval(header[0]) # number of records in an echo
  41. self.NE = eval(header[1]) # number of echoes
  42. self.DT = eval(header[2]) # in s
  43. self.TAUE = eval(header[3]) # in s
  44. self.RL = self.DT * self.NR # in s
  45. # pandas is about ~16 times faster than numpy for import (7 sec vs. .5 sec on tested example)
  46. # Parse in datafile, TODO update for newer format and remove hardcoded 3 in pandas
  47. #DATA = np.loadtxt(filename, comments="#", skiprows=ih)
  48. DATA = read_csv(filename, comment="#", sep='\t', skiprows=3, engine='c', header=None, names=('time','inphase','quadrature')).as_matrix()
  49. self.NS = int(np.shape(DATA)[0] / (self.NE * self.NR))
  50. # reshape DATA
  51. #self.DATA = np.reshape( DATA[:,1] + 1j*DATA[:,2], ( self.NS, self.NE, self.NR) )
  52. self.DATA = np.reshape( self.A2D(DATA[:,1]) + 1j*self.A2D(DATA[:,2]), ( self.NS, self.NE, self.NR) )
  53. self.WTIME = np.reshape( DATA[:,0], ( self.NS, self.NE, self.NR) )
  54. self.TIMES = np.arange(self.DT, self.RL + 1e-8, self.DT) # 1 ms record, 1 us dwell
  55. self.emit(SIGNAL("plotRAW()"))
  56. self.emit(SIGNAL("enableDSP()"))
  57. self.emit(SIGNAL("doneStatus()"))
  58. def evaluate(self, value):
  59. if len(value) == 1:
  60. if value[0][-1].isdigit():
  61. return eval(value[0])
  62. else:
  63. if value[0][-1] == "u": # micro
  64. return eval(value[0][0:-1]) * 1e-6
  65. elif value[0][-1] == "M": # mega
  66. return eval(value[0][0:-1]) * 1e6
  67. elif value[0][-1] == "s": # second
  68. return eval(value[0][0:-1])
  69. else:
  70. print("DOOM!")
  71. exit(1)
  72. else:
  73. # TODO
  74. # logic here is a little lacking, but current tool does not do anything
  75. # sneaky with array values, could break in the future though
  76. vals = []
  77. for val in value:
  78. vals.append( eval(val) )
  79. return vals
  80. def loadORSFile(self, fname):
  81. """ Loads ORS files with extensive header info
  82. """
  83. parameters = False
  84. data = False
  85. pars = {}
  86. with open(fname) as f:
  87. iline = 0
  88. for line in f:
  89. if len(line.split()) == 0:
  90. pass # empty do nothing
  91. else:
  92. if line.split()[0] == "[ps]":
  93. parameters = False
  94. # these entries specify the chip programming
  95. if parameters:
  96. # proc header info
  97. split,comment = line.split("#")
  98. var, val = split.split("=")
  99. val = val.split(";")[0].split()
  100. pars[var.strip()] = self.evaluate(val)
  101. if line.split()[0] == "[parameters]":
  102. parameters = True
  103. if data:
  104. if line[0] == "#":
  105. pass
  106. else:
  107. header = line.split()
  108. iline += 1
  109. break
  110. if line.split()[0] == "[data]":
  111. data = True
  112. iline += 1
  113. self.NR = eval(header[0]) # number of records in an echo
  114. self.NE = eval(header[1]) # number of echoes
  115. self.DT = eval(header[2]) # in s
  116. self.TAUE = eval(header[3]) # in s
  117. self.RL = self.DT * self.NR # in s
  118. # Pandas is much faster than numpy
  119. DATA = read_csv(fname, comment="#", sep='\t', skiprows=iline, engine='c', header=None, names=('time','inphase','quadrature')).as_matrix()
  120. self.NS = int(np.shape(DATA)[0] / (self.NE * self.NR))
  121. # reshape DATA
  122. #self.DATA = np.reshape( DATA[:,1] + 1j*DATA[:,2], ( self.NS, self.NE, self.NR) )
  123. self.DATA = np.reshape( self.A2D(DATA[:,1], pars) + 1j*self.A2D(DATA[:,2], pars), ( self.NS, self.NE, self.NR) )
  124. self.WTIME = np.reshape( DATA[:,0], ( self.NS, self.NE, self.NR) )
  125. self.TIMES = np.arange(self.DT, self.RL + 1e-8, self.DT) # 1 ms record, 1 us dwell
  126. self.emit(SIGNAL("plotRAW()"))
  127. self.emit(SIGNAL("enableDSP()"))
  128. self.emit(SIGNAL("doneStatus()"))
  129. # analogue to digital conversion, 'don't quote me on this -- kevin'
  130. def A2D(self, i, pars):
  131. """ Actual voltage is not known. This routine will need to be expaned to account for
  132. calibration datasets. Or perhaps the inversion and fitting routines instead.
  133. TODO, if prestacked data, we need to account for this.
  134. """
  135. #return (i/(.5 * (2.**16)))
  136. #cal = 7.85 # TODO query for this
  137. cal = 5.51 # 1/2 full
  138. return i / ( pars['n'] * pars['cal'] )
  139. def loadJavelinLog(self, filename):
  140. """Loads processed Javelin data"""
  141. #print("Javelin time")
  142. import scipy.io as sio
  143. Data = sio.loadmat(filename)
  144. self.T2T = Data['time'][:,0]
  145. self.TAUE = self.T2T[3] - self.T2T[2]
  146. # Better sigma estimate needed
  147. self.sigma = 2.5
  148. self.T2N = np.random.normal( 0, self.sigma, len(self.T2T) )
  149. self.NS = 1
  150. self.T2D = self.T2N + 1j*Data['se_vector_wc'][:,56] # depth 0
  151. # What about each depth? How is modelling handled.
  152. def loadVCMFile(self, filename):
  153. #print("LOADING VCM file")
  154. Data = np.loadtxt(filename, comments = "#")
  155. self.T2T = Data[:,0]
  156. self.sigma = eval(open( filename, 'r' ).readline().split()[1])
  157. self.T2N = np.random.normal( 0, self.sigma, len(self.T2T) )
  158. self.T2D = self.T2N + 1j*Data[:,1]
  159. self.NS = 1
  160. self.TAUE = self.T2T[3] - self.T2T[2]
  161. #self.emit(SIGNAL("plotLOG10ENV()"))
  162. #self.emit(SIGNAL("enableDSP()"))
  163. #self.emit(SIGNAL("enableINV()"))
  164. #self.emit(SIGNAL("doneStatus()"))
  165. def loadTextFile(self, filename, NS, NE, NR):
  166. """ Loads older format files without header info """
  167. print ("OLD TEXT FORMAT DETECTED, not SUPPORTED RIGHT NOW")
  168. exit()
  169. pass
  170. def batchThread(self, pars, tag, bload=True, bwindow=True, benvelope=True, bphase=True, bgate=True):
  171. #print("batdacheing", load, window, envelope, phase, gate)
  172. fileName, fileExtension = os.path.splitext( str(tag) )
  173. if bload:
  174. #print ("loading", tag)
  175. if fileExtension == ".ors":
  176. self.loadORSFile(str(tag))
  177. elif fileExtension == ".vcm":
  178. self.loadVCMFile(str(tag))
  179. if bwindow:
  180. #print("windowing", pars['ftype'], pars['winWidth'])
  181. self.WindowAndStack( pars['ftype'] , pars['winTrimLeft'], pars['winTrimRight'] )
  182. if benvelope:
  183. #print("enveloping", pars['offset'])
  184. self.T2EnvelopeDetect( pars['offset'] )
  185. if bphase:
  186. self.CorrectedAmplitude( 0,0 )
  187. if bgate:
  188. #print("gating", pars['gpd'])
  189. self.gateIntegrate( pars['gpd'], pars['stackEfficiency'] )
  190. self.emit(SIGNAL("doneStatus()"))
  191. def WindowAndStack(self, ftype="Blackman-Harris", triml=0, trimr=0):
  192. # TODO, consider lowpass filter and interpolate to higher DT in order to minimize window amplitude effects
  193. # what would cost of such an operation be?
  194. # TODO, should record be overwritten? I don't think so. I think we should make a deep copy here and then you can always fall back.
  195. # Trim ends to centre echo, this many are REMOVED
  196. #triml = 18
  197. #trimr = 15
  198. self.TIMES = self.TIMES[triml:-(trimr)]
  199. self.NR -= triml + trimr
  200. if ftype == "Blackman-Harris":
  201. window = signal.blackmanharris(self.NR)
  202. elif ftype == "Blackman":
  203. window = signal.blackman(self.NR)
  204. elif ftype == "Hann":
  205. window = signal.hann(self.NR)
  206. elif ftype == "sinc":
  207. window = np.ones( (self.NR) )
  208. self.DATAP = self.DATA[:,:,triml:-trimr].copy()
  209. for istk in range(0,self.NS):
  210. for ie in range(0,self.NE,1):
  211. self.DATAP[istk, ie,:] *= window # signal.blackmanharris(self.NR)
  212. #pass
  213. #self.DATA[istk, ie, :] *= signal.hann(self.NR)
  214. #self.DATA[istk, ie, :] *= signal.blackman(self.NR)
  215. #DATA2[istk, ie, :] *= signal.blackmanharris(self.NR)
  216. self.emit(SIGNAL("updateProgress(int)"), (int)(1e2*istk/self.NS))
  217. self.DATASTACK = np.average(self.DATAP, axis=0)
  218. if (self.NS > 1):
  219. # for noise calculation, destructively average phase-cycled echoes
  220. self.DATAP[1::2,:,:] *= -1.
  221. self.NOISE = np.average(self.DATAP, axis=0)
  222. self.NOISE -= np.mean(self.NOISE) # The noise can be biased? Don't understand why!!
  223. #self.DATA[0::2,:,:] *= -1.
  224. self.emit(SIGNAL("plotWIN()"))
  225. self.emit(SIGNAL("enableDSP()"))
  226. self.emit(SIGNAL("doneStatus()"))
  227. def ResampleEnvelopeData(self, Q, Mask=0, GatesPD=0):
  228. """ Mirror the data and then apply window function and resample in Fourier domain at twice the number.
  229. Then take half. Should be better.
  230. GatesPD is Gates Per Decade
  231. """
  232. tt = np.concatenate( [self.T2T[Mask::], self.T2T[Mask::][-1]+self.T2T[Mask::] ] )
  233. xx = np.concatenate( [self.T2D[Mask::][::-1] , self.T2D[Mask::] ] )
  234. xx,tx = signal.resample( xx, 2*Q, t=tt)#, window='hanning')
  235. self.T2D = xx[0:Q][::-1]
  236. #self.T2D = self.T2D[::-1]
  237. xx = np.concatenate( [self.T2N[::-1] , self.T2N] )
  238. xx,tx = signal.resample( xx, 2*Q, t=tt)#, window='hanning')
  239. self.T2N = xx[0:Q][::-1]
  240. # don't reverse times dummy!
  241. self.T2T = tx[0:Q] #[::-1] ## <-- No! bad Trevor
  242. self.NE = len(self.T2D)
  243. def computeSTD(self):
  244. if self.NS >= 2:
  245. self.sigma = np.std( np.imag(self.T2N) )
  246. else:
  247. s0 = np.std( np.real(self.T2D) )
  248. self.sigma = s0
  249. def computeSTDBurst(self):
  250. if self.NS >= 2:
  251. self.sigmaBurst = np.std( np.imag(self.T2Nb) )
  252. else:
  253. N = np.real(self.T2Db)
  254. s0 = np.std( np.real(self.T2Db) )
  255. sr = ( np.std( (np.exp(-1j*self.zeta)*N).real) )
  256. si = ( np.std( (np.exp(-1j*self.zeta)*N).imag) )
  257. self.sigmaBurst = s0 #+ (sr+si-s0)/2.
  258. #self.T2Db.real *= self.sigmaBurst/s0 #nr #(nr+ni)/2.
  259. #self.sigmaBurst = np.std( np.real(self.T2Db) )
  260. def gateIntegrateBurst(self, gpd, stackEfficiency):
  261. self.computeSTDBurst()
  262. # calculate total number of decades
  263. nd = np.log10(self.T2Tb[-1]/self.T2Tb[0]) #np.log10(self.T2T[-1]) - np.log10(self.T2T[-1])
  264. tdd = np.logspace( np.log10(self.T2Tb[0]), np.log10(self.T2Tb[-1]), (int)(gpd*nd)+1, base=10, endpoint=True)
  265. tdl = tdd[0:-1] # these are left edges
  266. tdr = tdd[1::] # these are left edges
  267. td = (tdl+tdr) / 2. # window centres
  268. htd = np.zeros( len(td), dtype=complex )
  269. htn = np.zeros( len(td), dtype=complex )
  270. isum = np.zeros( len(td), dtype=int )
  271. Vars = np.ones( len(td) ) * self.sigmaBurst**2 #* .15**2
  272. ii = 0
  273. for itd in range(len(self.T2Tb)):
  274. if ( self.T2Tb[itd] > tdr[ii] ):
  275. if (ii < len(td)-1):
  276. ii += 1
  277. else:
  278. break
  279. isum[ii] += 1
  280. htd[ii] += self.T2Db[ itd ]
  281. #htn[ii] += self.T2N[ itd ]
  282. Vars[ii] += self.sigmaBurst**2
  283. #self.emit(SIGNAL("updateProgress(int)"), (int)(1e2*itd/len(self.T2T)))
  284. # The 1.5 should be 2 in the theory, the 1.5 compensates for correlated noise, etc. Bascically
  285. # we don't see noise reduce like ideal averaging.
  286. self.sigmaBurst = np.sqrt( Vars * (1/isum)**stackEfficiency ) # Var (aX) = a^2 Var(x)
  287. # Bootstrap
  288. bootstrap = False
  289. if bootstrap:
  290. Means, isumbs = self.bootstrapWindows(np.real(self.T2Db), 300, isum)
  291. ts,sm = self.smooth( isumbs, np.std(Means, axis=1) )
  292. for it in range(ii):
  293. self.sigmaBurst[it] = sm[isum[it]]
  294. # RESET times where isum == 1
  295. ii = 0
  296. while (isum[ii] == 1):
  297. td[ii] = self.T2Tb[ii]
  298. ii += 1
  299. htd /= isum
  300. htn /= isum
  301. self.T2Tb = td
  302. self.T2Db = htd
  303. def gateIntegrate(self, gpd, stackEfficiency):
  304. """
  305. Gate integrate the signal to gpd, gates per decade
  306. """
  307. self.computeSTD()
  308. self.gpd = gpd
  309. if self.burst:
  310. self.gateIntegrateBurst(gpd, stackEfficiency)
  311. print ("BBBBBUUURSSSSTTTTTTTTTTTT DOOOM FIXME MRProc.py!!!!!!!!!!!")
  312. # I need to apply gate shifting logic to burst data as well
  313. #########################################################################
  314. # use artificial time gates so that early times are fully captured at gpd
  315. # very minor correction for borehole data.
  316. T2T0 = self.T2T[0]
  317. T2TD = self.T2T[0] - (self.T2T[1]-self.T2T[0])
  318. self.T2T -= T2TD
  319. # calculate total number of decades
  320. nd = np.log10(self.T2T[-1]/self.T2T[0]) #np.log10(self.T2T[-1]) - np.log10(self.T2T[-1])
  321. tdd = np.logspace( np.log10(self.T2T[0]), np.log10(self.T2T[-1]), (int)(gpd*nd)+1, base=10, endpoint=True)
  322. tdl = tdd[0:-1] # these are left edges
  323. tdr = tdd[1::] # approximate left edges
  324. td = (tdl+tdr) / 2. # approximate window centres
  325. htd = np.zeros( len(td), dtype=complex )
  326. htn = np.zeros( len(td), dtype=complex )
  327. isum = np.zeros( len(td), dtype=int )
  328. Vars = np.zeros( len(td) )
  329. ii = 0
  330. for itd in range(len(self.T2T)):
  331. if ( self.T2T[itd] > tdr[ii] ):
  332. if (ii < len(td)-1):
  333. ii += 1
  334. # correct window edges to centre about data
  335. tdr[ii-1] = (self.T2T[itd-1]+self.T2T[itd])*.5
  336. tdl[ii ] = (self.T2T[itd-1]+self.T2T[itd])*.5
  337. else:
  338. break
  339. isum[ii] += 1
  340. htd[ii] += self.T2D[ itd ]
  341. #htn[ii] += self.T2N[ itd ]
  342. Vars[ii] += self.sigma**2
  343. self.emit(SIGNAL("updateProgress(int)"), (int)(1e2*itd/len(self.T2T)))
  344. # Correct window centres
  345. td = (tdl+tdr) / 2. # actual window centres
  346. # Theoretical gates
  347. # The 1.5 should be 2 in the theory, the 1.5 compensates for correlated noise, etc. Bascically
  348. # we don't see noise reduce like ideal averaging.
  349. self.sigma = np.sqrt(Vars * (1./isum)**stackEfficiency ) # Var (aX) = a^2 Var(x)
  350. # Bootstrap
  351. bootstrap = False
  352. #bs = open("bootstrap.dat", "a")
  353. if bootstrap:
  354. print("Bootstrappin")
  355. nboot=20000
  356. Means, isumbs = self.bootstrapWindows(np.real(self.T2D), nboot, isum[isum!=1], adapt=True)
  357. #Means, isumbs = self.bootstrapWindows(np.real(self.T2D), nboot, np.arange(1,isum[-1]+1), adapt=True)
  358. # MAD measure
  359. c = stats.norm.ppf(3/4.)
  360. mad = np.ma.median(np.ma.abs(Means), axis=1)/c
  361. #nboot2 = nboot//isumbs
  362. #std = np.ma.std(Means, axis=1, ddof=2) # unbiased
  363. #iqr = stats.iqr( Means, axis=1 )
  364. # Smoothing spline
  365. #ts,sm = self.smooth(isumbs, std)
  366. #ts,sm = self.smooth(isumbs, std)
  367. #for ii, s in enumerate(mad):
  368. # bs.write( str( (s-self.sigma[isum!=1][ii]) / self.sigma[isum!=1][ii] ) + " " )
  369. #for item in (sm[isum[isum!=1]-1]-self.sigma[isum!=1])/self.sigma[isum!=1]:
  370. # bs.write( str(item) + " ")
  371. #bs.write("\n")
  372. #ts,sm = self.smooth( isumbs, np.std(Means, axis=1))
  373. #bias = np.average(Means,axis=1)
  374. ii = 0
  375. for it in range(1, len(self.sigma)):
  376. #self.sigma[it] = sm[isum[it]-1]
  377. if isum[it] != 1:
  378. #print ( "assigning mad", it , mad[ii] )
  379. self.sigma[it] = mad[ii]
  380. ii += 1
  381. # RESET times where isum == 1
  382. ii = 0
  383. while (isum[ii] == 1):
  384. td[ii] = self.T2T[ii]
  385. ii += 1
  386. htd /= isum
  387. htn /= isum
  388. self.T2T = td + T2TD
  389. self.T2D = htd
  390. self.emit(SIGNAL("plotLOG10ENV()"))
  391. self.emit(SIGNAL("enableDSP()"))
  392. self.emit(SIGNAL("enableINV()"))
  393. self.emit(SIGNAL("doneStatus()"))
  394. def bootstrapWindows(self, N, nboot, isum, adapt=False):
  395. """ Bootstraps noise as a function of gate width
  396. """
  397. nc = np.shape(N)[0]
  398. Means = {}
  399. if adapt:
  400. Means = -9999*np.ones((len(isum), nboot//isum[0]))
  401. for ii, nwin in enumerate(isum):
  402. for iboot in range(nboot//isum[ii]):
  403. cs = np.random.randint(0,nc-nwin)
  404. Means[ii,iboot] = np.mean( N[cs:cs+nwin] )
  405. #Means[ii,iboot] = np.mean( np.random.normal(0, self.sigma[0], nwin))
  406. Means = np.ma.masked_less(Means, -9995)
  407. else:
  408. Means = np.zeros((len(isum), nboot))
  409. for ii, nwin in enumerate(isum):
  410. for iboot in range(nboot):
  411. cs = np.random.randint(0,nc-nwin)
  412. Means[ii,iboot] = np.mean( N[cs:cs+nwin] )
  413. #Means[ii,iboot] = np.mean( np.random.normal(0, self.sigma[0], nwin))
  414. return Means, np.array(isum)
  415. def smooth(self, x, y):
  416. w = np.ones(len(x))
  417. #w[0] = 1. # force first time gate
  418. xs = np.arange(1, x[-1]+1, .1) # resample
  419. #s = UnivariateSpline( np.log(x), np.log(y), s=2.5, w=w)
  420. s = UnivariateSpline( np.log(x), np.log(y), s=4.5, w=w)
  421. #s = UnivariateSpline( x, y, s=2.5, w=w)
  422. ys = s(np.log(xs))
  423. #ys = s(xs)
  424. # resample
  425. ys = ys[0::10]
  426. xs = xs[0::10]
  427. return xs,np.exp(ys)
  428. def T2EnvelopeDetect(self, offset=0):
  429. self.NRS = (int)(2**12) # Number of zeros to pad with, total record length
  430. # Zero pad
  431. DATASTACKZPF = np.append( self.DATASTACK, np.zeros( (self.NE, self.NRS - self.NR) ), axis=1 )
  432. if (self.NS > 1):
  433. NOISEZPF = np.append( self.NOISE, np.zeros( (self.NE, self.NRS - self.NR) ), axis=1 )
  434. # TODO check for proper split
  435. SPLIT = self.NR//2 + offset # 20 # 15 #2
  436. # containers for array
  437. self.T2D = np.zeros( (self.NE), dtype='complex' )
  438. self.T2N = np.zeros( (self.NE), dtype='complex' )
  439. self.T2T = np.zeros( (self.NE) )
  440. # frequency-domain processing
  441. for ie in range(0,self.NE):
  442. COLOUR = np.array([ie/self.NE, 0., 1.-ie/self.NE]) # for plots
  443. ###############################################################
  444. # fold data need to reverse some of these
  445. #FOLD = ((DATASTACK[ie,0:NR/2] + DATASTACK[ie,:][::-1][0:NR/2])/2.)[::-1]
  446. #plt.figure(223)
  447. #plt.plot( self.TAUE*(float)(ie) + self.TIMES, np.imag(self.DATASTACK[ie]), color=COLOUR)
  448. #plt.plot( self.TAUE*(float)(ie) + self.TIMES, np.imag(self.DATA[0,ie]), color='green')
  449. #plt.plot( self.TAUE*(float)(ie) + self.TIMES, np.imag(self.DATA[1,ie]), color='black')
  450. #plt.plot( self.TAUE*(float)(ie) + self.TIMES, np.imag(self.DATA[2,ie]), color='cyan')
  451. #plt.plot( self.TAUE*(float)(ie) + self.TIMES, np.imag(self.DATA[3,ie]), color='purple')
  452. #plt.plot( TAUE*(1e-6)*(float)(ie) + TIMES[SPLIT], np.imag(DATASTACK[ie])[SPLIT], 'o', markersize=6, color=COLOUR)
  453. #plt.plot( TAUE*DT*ie + TIMES, np.real(DATASTACK[ie]), color=COLOUR[::-1])
  454. ##################################################################
  455. # do frequency domain workflow
  456. # fold to DC
  457. ND = np.append( DATASTACKZPF[ie,:][SPLIT::], DATASTACKZPF[ie,:][0:SPLIT] )
  458. #X = 1e-2*((self.RL/2.)/self.NR) * np.fft.fft(ND, n=self.NRS) # Why 1e-2?
  459. X = 1e4*(( (self.DT*self.NRS) /2.)/self.NRS) * np.fft.fft(ND, n=self.NRS) # Why 1e-2?
  460. #X = (self.DT) * np.fft.fft(ND, n=self.NRS) # Why 1e-2?
  461. if (self.NS > 1):
  462. NND = np.append( NOISEZPF[ie,:][SPLIT::], NOISEZPF[ie,:][0:SPLIT] )
  463. #XN = 1e-2*((self.RL/2.)/self.NR) * np.fft.fft(NND, n=self.NRS) # Why 1e-2?
  464. XN = 1e4*(( (self.DT*self.NRS) /2.)/self.NRS) * np.fft.fft(NND, n=self.NRS) # Why 1e-2?
  465. self.T2N[ie] = XN[0] #complex( np.real(X[0]), np.imag(X[0]) )
  466. ##################################################################
  467. # Save T2 records
  468. self.T2D[ie] = X[0] #complex( np.real(X[0]), np.imag(X[0]) )
  469. self.T2T[ie] = (1.+ie)*(self.TAUE)
  470. self.emit(SIGNAL("updateProgress(int)"), (int)(1e2*ie/self.NE))
  471. self.computeSTD()
  472. self.emit(SIGNAL("plotENV()"))
  473. self.emit(SIGNAL("enableDSP()"))
  474. self.emit(SIGNAL("enableINV()"))
  475. self.emit(SIGNAL("doneStatus()"))
  476. def CorrectedAmplitude(self, joe=0, frank=0):
  477. """ Phase corrects the record, joe and Frank are empty (ignored) and are workarounds for threading.
  478. """
  479. #[E0,df,phi,T2] = quadratureDetect(Q, I, tt)
  480. # mask first few echoes, just for detection
  481. NM = 2
  482. #[conv,E0,df,phi,T2] = quadratureDetect(np.imag(self.T2D[NM::]), np.real(self.T2D[NM::]), self.T2T[NM::], False, False, False)
  483. #[conv,E0,df,phi,T2] = quadratureDetect(np.imag(self.T2D[NM::]), np.real(self.T2D[NM::]), self.T2T[NM::], False, False, True)
  484. # #CorrectFreq=False, BiExp=False, CorrectDC=False)
  485. [conv,E0,df,phi,T2] = quadratureDetect2(np.imag(self.T2D[NM::]), np.real(self.T2D[NM::]), self.T2T[NM::])
  486. print("df", df, "phi", phi)
  487. self.zeta1 = phi
  488. #if conv: # failed convergence
  489. # [conv,E0,df,phi,T2] = quadratureDetect(np.imag(self.T2D[NM::]), np.real(self.T2D[NM::]), self.T2T[NM::], False, False, False)
  490. self.T2D = self.RotateAmplitude(np.real(self.T2D), np.imag(self.T2D), phi, df, self.T2T)
  491. self.computeSTD()
  492. #print(str(conv) + "\t" + str(phi)) #, file = "phase.dat")
  493. if self.burst:
  494. #[conv,E0,df,phi,T2] = quadratureDetect(np.imag(self.T2Db), np.real(self.T2Db), self.T2Tb, False, False, False)
  495. [conv,E0,df,phi,T2] = quadratureDetect2(np.imag(self.T2Db), np.real(self.T2Db), self.T2Tb) #, False, False, False)
  496. #if conv: # failed convergence
  497. # [conv,E0,df,phi,T2] = quadratureDetect(np.imag(self.T2Db), np.real(self.T2Db), self.T2Tb, False, False, False)
  498. env = E0*np.exp(-np.array(self.T2Tb)/T2)
  499. self.T2Db = self.RotateAmplitude(np.real(self.T2Db), np.imag(self.T2Db), phi, df, self.T2Tb)
  500. self.computeSTDBurst()
  501. self.emit(SIGNAL("plotENV()"))
  502. self.emit(SIGNAL("enableDSP()"))
  503. self.emit(SIGNAL("enableINV()"))
  504. self.emit(SIGNAL("doneStatus()"))
  505. def RotateAmplitude(self, X, Y, zeta, df, t):
  506. self.zeta = zeta
  507. self.df = df
  508. V = X + 1j*Y
  509. return np.abs(V) * np.exp(1j * (np.angle(V) - zeta - 2.*np.pi*df*t))
  510. def MonoFit(self, mask=2, intercept="False"):
  511. component='imag'
  512. #print "MonoRes", component, mask, intercept
  513. a0 = 0
  514. if intercept=="True":
  515. if component == "imag":
  516. [a0,b0,rt20] = regressCurve(np.imag(self.T2D[mask::]), self.T2T[mask::], intercept=True)
  517. if component == "real":
  518. [a0,b0,rt20] = regressCurve(np.real(self.T2D[mask::]), self.T2T[mask::], intercept=True)
  519. rfit = r"$\tilde{T}_2$= %4d [ms]"%(1e3*rt20)
  520. else:
  521. if component == "imag":
  522. [b0,rt20] = regressCurve(np.imag(self.T2D[mask::]), self.T2T[mask::], intercept=False) #, intercept=True) #
  523. if component == "real":
  524. [b0,rt20] = regressCurve(np.real(self.T2D[mask::]), self.T2T[mask::], intercept=False) #, intercept=True) #
  525. rfit = r"$\tilde{T}_2$= %4d [ms]"%(1e3*rt20)
  526. extrapolate = False #False #True #False
  527. if extrapolate:
  528. Textr = 3
  529. nede = np.log10(Textr/self.T2T[-1]) #np.log10(self.T2T[-1]) - np.log10(self.T2T[-1])
  530. tex = np.logspace( np.log10(self.T2T[-1]), np.log10(Textr), (int)(self.gpd*nede)+1, base=10, endpoint=True)
  531. dex = 1j* (a0 + b0*np.exp(-np.array(tex)/rt20))
  532. self.T2T = np.concatenate( (self.T2T, tex) )
  533. self.T2D = np.concatenate( (self.T2D, dex) )
  534. self.sigma = np.concatenate( (self.sigma, self.sigma[-1]*np.ones(len(dex)) ) )
  535. #env = b0*np.exp(-np.array(self.T2T)/rt20)
  536. #self.MonoRes = { 'env':env, 'rfit':rfit, 'b0':b0, 'rt20' :rt20 } #[a0,b0,rt20] }
  537. env = a0 + b0*np.exp(-np.array(self.T2T)/rt20)
  538. self.MonoRes = { 'env':env, 'rfit':rfit, 'a0':a0, 'b0':b0, 'rt20' :rt20, 't':self.T2T } #[a0,b0,rt20] }
  539. self.emit(SIGNAL("plotMONO()"))
  540. self.emit(SIGNAL("enableDSP()"))
  541. self.emit(SIGNAL("enableINV()"))
  542. self.emit(SIGNAL("doneStatus()"))
  543. def BiFit(self, mask=2, intercept="False"):
  544. component='imag'
  545. #print "MonoRes", component, mask, intercept
  546. a0 = 0
  547. if intercept=="True":
  548. if component == "imag":
  549. [a0,b0,rt20,bb0,rrt20] = regressCurve2(np.imag(self.T2D[mask::]), self.T2T[mask::], sigma2=self.sigma[mask::], intercept=True)
  550. if component == "real":
  551. [a0,b0,rt20,bb0,rrt20] = regressCurve2(np.real(self.T2D[mask::]), self.T2T[mask::], sigma2=self.sigma[mask::], intercept=True)
  552. rfit = r"$\tilde{T}_2$= %4d [ms]"%(1e3*rt20)
  553. else:
  554. if component == "imag":
  555. [b0,rt20,bb0,rrt20] = regressCurve2(np.imag(self.T2D[mask::]), self.T2T[mask::], sigma2=self.sigma[mask::], intercept=False) #, intercept=True) #
  556. if component == "real":
  557. [b0,rt20,bb0,rrt20] = regressCurve2(np.real(self.T2D[mask::]), self.T2T[mask::], sigma2=self.sigma[mask::], intercept=False) #, intercept=True) #
  558. rfit = r"$\tilde{T}_2$= %4d [ms]"%(1e3*rt20)
  559. extrapolate = True # True #True
  560. if extrapolate:
  561. Textr = 3
  562. nede = np.log10(Textr/self.T2T[-1]) #np.log10(self.T2T[-1]) - np.log10(self.T2T[-1])
  563. tex = np.logspace( np.log10(self.T2T[-1]), np.log10(Textr), (int)(self.gpd*nede)+1, base=10, endpoint=True)
  564. dex = 1j* (a0 + b0*np.exp(-np.array(tex)/rt20) + bb0*np.exp(-np.array(tex)/rrt20) )
  565. self.T2T = np.concatenate( (self.T2T, tex) )
  566. self.T2D = np.concatenate( (self.T2D, dex) )
  567. self.sigma = np.concatenate( (self.sigma, self.sigma[-1]*np.ones(len(dex)) ) )
  568. #env = b0*np.exp(-np.array(self.T2T)/rt20)
  569. #self.MonoRes = { 'env':env, 'rfit':rfit, 'b0':b0, 'rt20' :rt20 } #[a0,b0,rt20] }
  570. env = a0 + b0*np.exp(-np.array(self.T2T)/rt20) + bb0*np.exp(-np.array(self.T2T)/rrt20)
  571. #print ("bi fits", a0, b0, rt20, bb0, rrt20)
  572. self.BiRes = { 'env':env, 'rfit':rfit, 'a0':a0, 'b0':b0, 'rt20' :rt20, 'bb0':bb0, 'rrt20':rrt20, 't':self.T2T } #[a0,b0,rt20] }
  573. self.emit(SIGNAL("plotBI()"))
  574. self.emit(SIGNAL("enableDSP()"))
  575. self.emit(SIGNAL("enableINV()"))
  576. self.emit(SIGNAL("doneStatus()"))
  577. def DistFit(self, dc=0, mask=0, nT2=30, lowT2=.005, hiT2=3.25, spacing="Log_10", Smooth="Smallest", betaScale=1):
  578. Time = pwcTime()
  579. Time.setT2(lowT2, hiT2, nT2, spacing)
  580. Time.setSampling( self.T2T[mask:] )
  581. Time.generateGenv()
  582. #self.burst = False
  583. #print ("MRProc.py ~~ DistFit self.burst override for publication ", self.burst)
  584. if self.burst:
  585. # fit burst data first
  586. Timeb = pwcTime()
  587. #Timeb.setT2(lowT2, hiT2, nT2, spacing)
  588. Timeb.T2Bins = np.copy(Time.T2Bins)
  589. #if len(Timeb.T2Bins) > 20:
  590. # Timeb.T2Bins = Timeb.T2Bins[0:20] # Subset
  591. Timeb.setSampling( self.T2Tb )
  592. Timeb.generateGenv()
  593. # invert burst data, always use smooth to encourage fast time inversion
  594. modelb = logBarrier(Timeb.Genv, np.imag(self.T2Db)-dc, Timeb.T2Bins, "lcurve", MAXITER=500, sigma=betaScale*self.sigmaBurst, alpha=1e6, smooth=Smooth)[0]
  595. modelb2 = np.zeros(nT2)
  596. modelb2[0:len(modelb)] += modelb
  597. #model = modelb2 # TODO remove
  598. # invert regular data independently, sum inversions
  599. #model = logBarrier(Time.Genv, np.imag(self.T2D[mask:])-dc, MAXITER=500, x_0=modelb, sigma=betaScale*self.sigma[mask:], alpha=1e10, smooth=Smooth)
  600. #model += modelb
  601. # Pass burt result as reference model
  602. model = logBarrier(Time.Genv, np.imag(self.T2D[mask:])-dc, Time.T2Bins, "lcurve", MAXITER=500, xr=modelb2, sigma=betaScale*self.sigma[mask:], alpha=1e6, smooth=Smooth)[0]
  603. else:
  604. model = logBarrier(Time.Genv, np.imag(self.T2D[mask:])-dc, Time.T2Bins, "lcurve", MAXITER=500, sigma=betaScale*self.sigma[mask:], alpha=1e6, smooth=Smooth)[0] #, smooth=True) #
  605. # forward model
  606. env = np.dot(Time.Genv, model)
  607. # NOTE this is not thread safe or even thread aware. Need to use MPI style message passing
  608. self.DistRes = { 'env':env, 'mod':model, 'Time':Time } #[a0,b0,rt20] } no dice
  609. def DistFitMP(self, q, tag, dc=0, mask=0, nT2=30, lowT2=.005, hiT2=3.25, spacing="Log_10", Smooth="Smallest", betaScale=1):
  610. """ Multi-processor version of DistFit using Queue. Note that Queue doesn't
  611. tolerate large messages. So instead we write them to file, using pickle!
  612. It's semi-absurd. But seems to function fast enough.
  613. """
  614. import pickle
  615. self.DistFit(dc, mask, nT2, lowT2, hiT2, spacing, Smooth, betaScale)
  616. self.DistRes['tag'] = tag
  617. pickle.dump( self.DistRes, open( str(tag)+".p", "wb" ) )
  618. q.put( tag )
  619. if __name__ == "__main__":
  620. fileName, fileExtension = os.path.splitext(sys.argv[1])
  621. if fileExtension == ".ors":
  622. ORS = MRProc()
  623. ORS.loadORSFile(sys.argv[1])
  624. ORS.WindowAndStack()
  625. ORS.T2EnvelopeDetect()
  626. #ORS.ResampleEnvelopeData(200, Mask=2, GatesPD=20) # resample size
  627. #print (len(ORS.T2D))
  628. #exit()
  629. #mono, rfit, [a0,b0,mt2] = ORS.MonoFit()
  630. mono, rfit, [b0,mt2] = ORS.MonoFit()
  631. #ca, caEnv, CAS = ORS.CorrectedAmplitude()
  632. mymask = 2
  633. env, mod, Time = ORS.DistFit(0, mymask) #a0)
  634. # Subtract DC term before 2D fit? Kind of kludgy but I don't see a way around it right now.
  635. plt.figure()
  636. #plt.plot(ORS.T2T, np.imag(ca), label='imag ca')
  637. #plt.plot(ORS.T2T, np.real(ca), label='real ca')
  638. #if ORS.NS < 2:
  639. # plt.plot(ORS.T2T, np.imag(ORS.T2D), label='imag')
  640. # plt.plot(ORS.T2T, np.real(ORS.T2D), label='real')
  641. #else:
  642. plt.errorbar(ORS.T2T, np.imag(ORS.T2D), fmt='o', yerr=ORS.sigma, label="data imag")
  643. plt.errorbar(ORS.T2T, np.real(ORS.T2D), fmt='o', yerr=np.std(ORS.T2D.real), label="data real")
  644. if ORS.NS > 2:
  645. plt.errorbar(ORS.T2T, np.imag(ORS.T2N), fmt='o', yerr=ORS.sigma, label="noise imag")
  646. plt.plot(ORS.T2T, mono, color='cyan', linewidth=3, label=rfit)
  647. #plt.plot(ORS.T2T, caEnv, label="car")
  648. plt.plot(ORS.T2T[mymask::], env, color='magenta', linewidth=3, label="dist")
  649. #plt.plot(ORS.T2T, a0+env)
  650. plt.legend()
  651. plt.savefig(sys.argv[1].strip(".ors")+"_caenv.pdf", facecolor='white', edgecolor='white', dpi=200)
  652. # Report RMS error of two results
  653. print("RMS error mono-exponential ", np.linalg.norm(np.imag(ORS.T2D) - mono) )
  654. print("RMS error dist-fit ", np.linalg.norm(np.imag(ORS.T2D[mymask::]) - env) )
  655. print("RMS error real channel ", np.linalg.norm(np.real(ORS.T2D) ) )
  656. if ORS.NS > 2:
  657. print("RMS error noise channel ", np.linalg.norm(np.imag(ORS.T2N) ) )
  658. plt.figure()
  659. plt.plot( ORS.T2T[mymask::], np.imag(ORS.T2D[mymask::]) - env )
  660. plt.title("residuals")
  661. #plt.show()
  662. #exit()
  663. plt.figure()
  664. plt.plot(Time.T2Bins, mod, color='purple', linewidth=2, label="recovered")
  665. #plt.plot(Time.T2Bins, guess, color='red', linewidth=2, label="guess")
  666. axmine = plt.gca()
  667. axmine.set_xlabel(r"$T_2$ [s]", color='black')
  668. axmine.set_ylabel(r"$A_0$ [rku]", color='black')
  669. axmine2 = plt.twinx()
  670. axmine2.set_ylabel(r"$A_0$ [rku]", color='black')
  671. theta = np.sum(mod)
  672. LogMeanT2 = np.exp(np.sum( mod * np.log(Time.T2Bins) ) / theta )
  673. #axmine2.plot(mt2, a0+b0, 'o', markersize=8, label="mono")
  674. axmine2.plot(mt2, b0, 'o', markersize=8, label="mono")
  675. #axmine2.plot(Time.T2Bins, guess, '-',color=colour[ic], linewidth=2, label="total")
  676. axmine2.plot(LogMeanT2, theta, 's', markersize=8, label="$T_{2ML}$")
  677. axmine2.set_ylim([0, axmine2.get_ylim()[1]])
  678. leg = plt.legend()
  679. plt.savefig(sys.argv[1].strip(".ors")+"_2DT2.pdf", facecolor='white', edgecolor='white', dpi=200)
  680. plt.figure(23)
  681. stats.probplot( np.imag(ORS.T2D) - mono , dist="norm", plot=plt) #, label="mono resid")
  682. stats.probplot( np.imag(ORS.T2D[mymask::]) - env , dist="norm", plot=plt) #, label="dist resid")
  683. if ORS.NS > 2:
  684. stats.probplot( np.imag(ORS.T2N), dist="norm", plot=plt) #, label="noise" )
  685. plt.savefig(sys.argv[1].strip(".ors")+"_qq.pdf")
  686. plt.show()
  687. exit()