Surface NMR processing and inversion GUI
Ви не можете вибрати більше 25 тем Теми мають розпочинатися з літери або цифри, можуть містити дефіси (-) і не повинні перевищувати 35 символів.

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675
  1. from akvo.tressel.SlidesPlot import *
  2. import numpy as np
  3. import sys
  4. import matplotlib.pyplot as plt
  5. import cmocean
  6. from pylab import meshgrid
  7. from akvo.tressel.logbarrier import *
  8. import yaml,os
  9. import multiprocessing
  10. import itertools
  11. from scipy.linalg import svd
  12. from matplotlib.backends.backend_pdf import PdfPages
  13. from matplotlib.colors import LogNorm
  14. from matplotlib.colors import LightSource
  15. from matplotlib.ticker import ScalarFormatter
  16. from matplotlib.ticker import MaxNLocator
  17. from matplotlib.ticker import AutoMinorLocator
  18. from matplotlib.ticker import LogLocator
  19. from matplotlib.ticker import FormatStrFormatter
  20. from matplotlib.colors import Normalize
  21. import cmocean
  22. from akvo.tressel.lemma_yaml import *
  23. from akvo.tressel import nonlinearinv as nl
  24. import pandas as pd
  25. import matplotlib.colors as colors
  26. # From https://stackoverflow.com/questions/18926031/how-to-extract-a-subset-of-a-colormap-as-a-new-colormap-in-matplotlib
  27. def truncate_colormap(cmap, minval=0.0, maxval=1.0, n=100):
  28. new_cmap = colors.LinearSegmentedColormap.from_list(
  29. 'trunc({n},{a:.2f},{b:.2f})'.format(n=cmap.name, a=minval, b=maxval),
  30. cmap(np.linspace(minval, maxval, n)))
  31. return new_cmap
  32. def buildKQT(K0,tg,T2Bins):
  33. """
  34. Constructs a QT inversion kernel from an initial amplitude one.
  35. """
  36. nlay, nq = np.shape(K0)
  37. nt2 = len(T2Bins)
  38. nt = len(tg)
  39. KQT = np.zeros( ( nq*nt,nt2*nlay), dtype=np.complex128 )
  40. for iq in range(nq):
  41. for it in range(nt):
  42. for ilay in range(nlay):
  43. for it2 in range(nt2):
  44. #KQT[iq*nt + it,ilay*nt2+it2] = K0[ilay,iq]*np.exp(-((10+tg[it])*1e-3)/(1e-3*T2Bins[it2]))
  45. KQT[iq*nt + it,ilay*nt2+it2] = K0[ilay,iq]*np.exp(-((10+tg[it])*1e-3)/(1e-3*T2Bins[it2]))
  46. return KQT
  47. def loadAkvoData(fnamein, chan):
  48. """ Loads data from an Akvo YAML file. The 0.02 is hard coded as the pulse length. This needs to be
  49. corrected in future kernel calculations. The current was reported but not the pulse length.
  50. """
  51. fname = (os.path.splitext(fnamein)[0])
  52. with open(fnamein, 'r') as stream:
  53. try:
  54. AKVO = (yaml.load(stream, Loader=yaml.Loader))
  55. except yaml.YAMLError as exc:
  56. print(exc)
  57. exit()
  58. Z = np.zeros( (AKVO.nPulseMoments, AKVO.Gated["Pulse 1"]["abscissa"].size ) )
  59. ZS = np.zeros( (AKVO.nPulseMoments, AKVO.Gated["Pulse 1"]["abscissa"].size ) )
  60. for q in range(AKVO.nPulseMoments):
  61. Z[q] = AKVO.Gated["Pulse 1"][chan]["Q-"+str(q) + " CA"].data
  62. if chan == "Chan. 1":
  63. ZS[q] = AKVO.Gated["Pulse 1"][chan]["STD"].data
  64. elif chan == "Chan. 2":
  65. ZS[q] = AKVO.Gated["Pulse 1"][chan]["STD"].data
  66. elif chan == "Chan. 3":
  67. ZS[q] = AKVO.Gated["Pulse 1"][chan]["STD"].data
  68. elif chan == "Chan. 4":
  69. ZS[q] = AKVO.Gated["Pulse 1"][chan]["STD"].data
  70. else:
  71. print("DOOM!!!")
  72. exit()
  73. #Z *= 1e-9
  74. #ZS *= 1e-9
  75. J = AKVO.Pulses["Pulse 1"]["current"].data
  76. J = np.append(J,J[-1]+(J[-1]-J[-2]))
  77. Q = AKVO.pulseLength[0]*J
  78. return Z, ZS, AKVO.Gated["Pulse 1"]["abscissa"].data, Q
  79. def catLayers(K0):
  80. K = np.zeros( (len(K0.keys()), len(K0["layer-0"].data)) , dtype=complex )
  81. for lay in range(len(K0.keys())):
  82. #print(K0["layer-"+str(lay)].data) # print (lay)
  83. K[lay] =K0["layer-"+str(lay)].data # print (lay)
  84. return 1e9*K # invert in nV
  85. def loadK0(fname):
  86. """ Loads in initial amplitude kernel
  87. """
  88. print("loading K0", fname)
  89. with open(fname) as f:
  90. K0 = yaml.load(f, Loader=yaml.Loader)
  91. K = catLayers(K0.K0)
  92. ifaces = np.array(K0.Interfaces.data)
  93. return ifaces, K
  94. #return ifaces, np.abs(K)
  95. def invertDelta(G, V_n, T2Bins, sig, alphastar):
  96. """ helper function that simply calls logBarrier, simplfies parallel execution
  97. """
  98. model = logBarrier(G, V_n, T2Bins, "Single", MAXITER=1, sigma=sig, alpha=alphastar, smooth="Smallest")
  99. return model
  100. def main():
  101. if (len (sys.argv) < 2):
  102. print ("akvoQT invertParameters.yaml")
  103. exit()
  104. with open(sys.argv[1], 'r') as stream:
  105. try:
  106. cont = (yaml.load(stream, Loader=yaml.Loader))
  107. except yaml.YAMLError as exc:
  108. print(exc)
  109. ###############################################
  110. # Load in data
  111. ###############################################
  112. V = []
  113. VS = []
  114. QQ = []
  115. tg = 0
  116. for dat in cont['data']:
  117. for ch in cont['data'][dat]['channels']:
  118. print("dat", dat, "ch", ch)
  119. v,vs,tg,Q = loadAkvoData(dat, ch)
  120. V.append(v)
  121. VS.append(vs)
  122. QQ.append(Q)
  123. for iv in range(1, len(V)):
  124. V[0] = np.concatenate( (V[0], V[iv]) )
  125. VS[0] = np.concatenate( (VS[0], VS[iv]) )
  126. V = V[0]
  127. VS = VS[0]
  128. ###############################################
  129. # Load in kernels
  130. ###############################################
  131. K0 = []
  132. for kern in cont["K0"]:
  133. ifaces,k0 = loadK0( kern )
  134. K0.append(k0)
  135. for ik in range(1, len(K0)):
  136. K0[0] = np.concatenate( (K0[0].T, K0[ik].T) ).T
  137. K0 = K0[0]
  138. #np.save("ifaces", ifaces)
  139. #exit()
  140. #plt.matshow(np.real(K0))
  141. #plt.show()
  142. #exit()
  143. ##############################################################
  144. # VERY Simple Sensitivity based calc. of noise per layer #
  145. # minimally useful, but retained for backwards compatibility #
  146. maxq = np.argmax(np.abs(K0), axis=1)
  147. maxK = .1 * np.abs(K0)[ np.arange(0,len(ifaces)-1), maxq ] # 10% water is arbitrary
  148. SNR = maxK / (VS[0][0])
  149. ###############################################
  150. # Build full kernel
  151. ###############################################
  152. T2Bins = np.logspace( np.log10(cont["T2Bins"]["low"]), np.log10(cont["T2Bins"]["high"]), cont["T2Bins"]["number"], endpoint=True, base=10)
  153. T2Bins2 = np.append( T2Bins, T2Bins[-1] + (T2Bins[-1]-T2Bins[-2]) )
  154. NT2 = len(T2Bins)
  155. KQT = np.real(buildKQT(np.abs(K0),tg,T2Bins))
  156. ###############################################
  157. # Linear Inversion
  158. ###############################################
  159. print("Calling inversion", flush=True)
  160. inv, ibreak, errn, phim, phid, mkappa, Wd, Wm, alphastar = logBarrier(KQT, np.ravel(V), T2Bins, "lcurve", MAXITER=150, sigma=np.ravel(VS), alpha=1e6, smooth="Smallest" )
  161. ################################
  162. # Summary plots, Data Space #
  163. ################################
  164. # TODO, need to clean this up for the case of multiple channels! Each channel should be a new row. It will be ugly, but important
  165. # TODO, loop over channels
  166. ich = 0
  167. for ch in cont['data'][dat]['channels']:
  168. figx = plt.figure( figsize=(pc2in(42.0),pc2in(22.)) )
  169. ax1 = figx.add_axes([.100, .15, .200, .70])
  170. ax2 = figx.add_axes([.325, .15, .200, .70]) # shifted to make room for shared colourbar
  171. axc1= figx.add_axes([.550, .15, .025, .70]) # shifted to make room for shared colourbar
  172. ax3 = figx.add_axes([.670, .15, .200, .70])
  173. axc2= figx.add_axes([.895, .15, .025, .70]) # shifted to make room for shared colourbar
  174. ax3.set_yscale('log')
  175. ax2.set_yscale('log')
  176. ax1.set_yscale('log')
  177. ax2.yaxis.set_ticklabels([])
  178. ax3.yaxis.set_ticklabels([])
  179. ax3.set_xscale('log')
  180. ax2.set_xscale('log')
  181. ax1.set_xscale('log')
  182. ax1.set_ylabel("Q (A $\cdot$ s)")
  183. ax1.set_xlabel("time (s)")
  184. ax2.set_xlabel("time (s)")
  185. ax3.set_xlabel("time (s)")
  186. #TT, QQQ = np.meshgrid(tg, np.ravel(QQ))
  187. TT, QQQ = np.meshgrid(tg, np.ravel(QQ[ich]))
  188. nq = np.shape(QQ[ich])[0] - 1 # to account for padding in pcolor
  189. nt = np.shape(tg)[0]
  190. ntq = nt*nq
  191. VV = V[ich*nq:ich*nq+nq,:] # slice this channel
  192. VVS = VS[ich*nq:ich*nq+nq,:] # slice this channel
  193. mmax = np.max(np.abs(VV))
  194. mmin = np.min(VV)
  195. obs = ax1.pcolor(TT, QQQ, VV, cmap=cmocean.cm.curl_r, vmin=-mmax, vmax=mmax, shading='auto') # pcolor edge not defined
  196. ax1.set_title("observed")
  197. pre = np.dot(KQT[ich*ntq:(ich+1)*ntq,:], inv)
  198. PRE = np.reshape( pre, np.shape(VV) )
  199. prem = ax2.pcolor(TT, QQQ, PRE, cmap=cmocean.cm.curl_r, vmin=-mmax, vmax=mmax,shading='auto' )
  200. ax2.set_title("predicted")
  201. cbar = plt.colorbar(prem, axc1)
  202. axc1.set_ylim( [np.min(VV), np.max(VV)] )
  203. cbar.outline.set_edgecolor(None)
  204. cbar.set_label('$V_N$ (nV)')
  205. DIFF = (PRE-VV) / VVS
  206. md = np.max(np.abs(DIFF))
  207. dim = ax3.pcolor(TT, QQQ, DIFF, cmap=cmocean.cm.balance, vmin=-md, vmax=md, shading='auto')
  208. ax3.set_title("misfit / $\widehat{\sigma}$")
  209. cbar2 = plt.colorbar(dim, axc2)
  210. #axc1.set_ylim( [np.min(V), np.max(V)] )
  211. cbar2.outline.set_edgecolor(None)
  212. cbar2.set_label('$V_N$ (nV)')
  213. #plt.colorbar(dim, ax3)
  214. figx.suptitle(ch + " linear Inversion")
  215. plt.savefig(ch + "dataspace.pdf")
  216. ich += 1
  217. ###############################################
  218. # Non-linear refinement!
  219. ###############################################
  220. nonLinearRefinement = cont['NonLinearRefinement']
  221. if nonLinearRefinement:
  222. KQTc = buildKQT(K0, tg, T2Bins)
  223. prec = np.abs(np.dot(KQTc, inv))
  224. phidc = np.linalg.norm(np.dot(Wd,prec-np.ravel(V)))**2
  225. print("PHID forward linear=", errn, "PHID forward nonlinear=", phidc/len(np.ravel(V)))
  226. res = nl.nonlinearinversion(inv, Wd, KQTc, np.ravel(V), Wm, alphastar )
  227. if res.success == True:
  228. INVc = np.reshape(res.x, (len(ifaces)-1,cont["T2Bins"]["number"]) )
  229. prec = np.abs(np.dot(KQTc, res.x))
  230. phidc = np.linalg.norm(np.dot(Wd,prec-np.ravel(V)))**2
  231. PREc = np.reshape( prec, np.shape(V) )
  232. print("PHID linear=", errn, "PHID nonlinear=", phidc/len(np.ravel(V)))
  233. while phidc/len(np.ravel(V)) > errn:
  234. phidc_old = phidc/len(np.ravel(V))
  235. #alphastar *= .9
  236. res = nl.nonlinearinversion(res.x, Wd, KQTc, np.ravel(V), Wm, alphastar )
  237. if res.success == True:
  238. INVc = np.reshape(res.x, (len(ifaces)-1,cont["T2Bins"]["number"]) )
  239. prec = np.abs(np.dot(KQTc, res.x))
  240. phidc = np.linalg.norm(np.dot(Wd,prec-np.ravel(V)))**2
  241. PREc = np.reshape( prec, np.shape(V) )
  242. print("PHID linear=", errn, "PHID nonlinear=", phidc/len(np.ravel(V)))
  243. else:
  244. break
  245. if phidc_old - phidc/len(np.ravel(V)) < 0.005:
  246. print("Not making progress reducing misfit in nonlinear refinement")
  247. break
  248. # Turn this into a nice figure w/ shared axes etc.
  249. # plt.matshow(PREc, cmap='Blues')
  250. # plt.gca().set_title("nonlinear predicted")
  251. # plt.colorbar()
  252. #
  253. # DIFFc = (PREc-V) / VS
  254. # md = np.max(np.abs(DIFF))
  255. # plt.matshow(DIFFc, cmap=cmocean.cm.balance, vmin=-md, vmax=md)
  256. # plt.gca().set_title("nonlinear misfit / $\widehat{\sigma}$")
  257. # plt.colorbar()
  258. ################################
  259. # Summary plots, Data Space #
  260. ################################
  261. ich = 0
  262. for ch in cont['data'][dat]['channels']:
  263. figx = plt.figure( figsize=(pc2in(42.0),pc2in(22.)) )
  264. ax1 = figx.add_axes([.100, .15, .200, .70])
  265. ax2 = figx.add_axes([.325, .15, .200, .70]) # shifted to make room for shared colourbar
  266. axc1= figx.add_axes([.550, .15, .025, .70]) # shifted to make room for shared colourbar
  267. ax3 = figx.add_axes([.670, .15, .200, .70])
  268. axc2= figx.add_axes([.895, .15, .025, .70]) # shifted to make room for shared colourbar
  269. ax3.set_yscale('log')
  270. ax2.set_yscale('log')
  271. ax1.set_yscale('log')
  272. ax2.yaxis.set_ticklabels([])
  273. ax3.yaxis.set_ticklabels([])
  274. ax3.set_xscale('log')
  275. ax2.set_xscale('log')
  276. ax1.set_xscale('log')
  277. ax1.set_ylabel("Q (A $\cdot$ s)")
  278. ax1.set_xlabel("time (s)")
  279. ax2.set_xlabel("time (s)")
  280. ax3.set_xlabel("time (s)")
  281. #TT, QQQ = np.meshgrid(tg, np.ravel(QQ))
  282. TT, QQQ = np.meshgrid(tg, np.ravel(QQ[ich]))
  283. nq = np.shape(QQ[ich])[0] - 1 # to account for padding in pcolor
  284. nt = np.shape(tg)[0]
  285. ntq = nt*nq
  286. VV = V[ich*nq:ich*nq+nq,:] # slice this channel
  287. VVS = VS[ich*nq:ich*nq+nq,:] # slice this channel
  288. mmax = np.max(np.abs(VV))
  289. mmin = np.min(VV)
  290. obs = ax1.pcolor(TT, QQQ, VV, cmap=cmocean.cm.curl_r, vmin=-mmax, vmax=mmax, shading='auto')
  291. ax1.set_title("observed")
  292. ## Here neds to change
  293. pre = np.abs(np.dot(KQTc[ich*ntq:(ich+1)*ntq,:], inv))
  294. PRE = np.reshape( pre, np.shape(VV) )
  295. prem = ax2.pcolor(TT, QQQ, PRE, cmap=cmocean.cm.curl_r, vmin=-mmax, vmax=mmax, shading='auto' )
  296. ax2.set_title("predicted")
  297. cbar = plt.colorbar(prem, axc1)
  298. axc1.set_ylim( [np.min(VV), np.max(VV)] )
  299. cbar.outline.set_edgecolor(None)
  300. cbar.set_label('$V_N$ (nV)')
  301. DIFF = (PRE-VV) / VVS
  302. md = np.max(np.abs(DIFF))
  303. dim = ax3.pcolor(TT, QQQ, DIFF, cmap=cmocean.cm.balance, vmin=-md, vmax=md, shading='auto')
  304. ax3.set_title("misfit / $\widehat{\sigma}$")
  305. cbar2 = plt.colorbar(dim, axc2)
  306. #axc1.set_ylim( [np.min(V), np.max(V)] )
  307. cbar2.outline.set_edgecolor(None)
  308. cbar2.set_label('$V_N$ (nV)')
  309. #plt.colorbar(dim, ax3)
  310. figx.suptitle(ch + " non-linear Inversion")
  311. plt.savefig(ch + "_NLdataspace.pdf")
  312. ich += 1
  313. ###############################################
  314. # Appraise DOI using simplified MRM
  315. ###############################################
  316. CalcDOI = cont['CalcDOI']
  317. if CalcDOI:
  318. pdf = PdfPages('resolution_analysis' + '.pdf' )
  319. MRM = np.zeros((len(ifaces)-1, len(ifaces)-1))
  320. # Build delta models
  321. DELTA = []
  322. for ilay in range(len(ifaces)-1):
  323. #for ilay in range(4):
  324. iDeltaT2 = len(T2Bins)//2
  325. deltaMod = np.zeros( (len(ifaces)-1, len(T2Bins)) )
  326. deltaMod[ilay][iDeltaT2] = 0.3
  327. dV = np.dot(KQT, np.ravel(deltaMod))
  328. #dinv, dibreak, derrn = logBarrier( KQT, dV, T2Bins, "single", MAXITER=1, sigma=np.ravel(VS), alpha=alphastar, smooth="Smallest" )
  329. #output = invertDelta(KQT, dV, T2Bins, np.ravel(VS), alphastar)
  330. DELTA.append(dV)
  331. print("Performing resolution analysis in parallel, printed output may not be inorder.", flush=True)
  332. with multiprocessing.Pool() as pool:
  333. invresults = pool.starmap(invertDelta, zip(itertools.repeat(KQT), DELTA, itertools.repeat(T2Bins), itertools.repeat(np.ravel(VS)), itertools.repeat(alphastar) ))
  334. # invresults = pool.starmap(logBarrier, zip(itertools.repeat(KQT), DELTA, itertools.repeat(T2Bins), itertools.repeat('single'), \
  335. # itertools.repeat('MAXITER=1'), itertools.repeat(np.ravel(VS)), itertools.repeat(alphastar))) #, itertools.repeat(u'smooth=\'Smallest\'')) )
  336. # This could be parallelized
  337. for ilay in range(len(ifaces)-1):
  338. # invert
  339. #dinv, dibreak, derrn = logBarrier(KQT, dV, T2Bins, "single", MAXITER=1, sigma=np.ravel(VS), alpha=alphastar, smooth="Smallest" )
  340. #print("Sum dinv from", str(ifaces[ilay]), "to", str(ifaces[ilay+1]), "=", np.sum(dinv))
  341. dinv, dibreak, derrn = invresults[ilay]
  342. DINV = np.reshape(dinv, (len(ifaces)-1,cont["T2Bins"]["number"]) )
  343. MRM[ilay,:] = np.sum(DINV, axis=1)
  344. Y,X = meshgrid( ifaces, T2Bins2 )
  345. fig = plt.figure( figsize=(pc2in(20.0),pc2in(22.)) )
  346. ax1 = fig.add_axes( [.2,.15,.6,.7] )
  347. im = ax1.pcolor(X, Y, DINV.T, cmap=cmocean.cm.tempo, shading='auto')
  348. ax1.plot( T2Bins[iDeltaT2], (ifaces[ilay]+ifaces[ilay+1])/2, 's', markersize=6, markeredgecolor='black') #, markerfacecolor=None )
  349. im.set_edgecolor('face')
  350. ax1.set_xlabel(u"$T_2^*$ (ms)")
  351. ax1.set_ylabel(u"depth (m)")
  352. ax1.set_xlim( T2Bins2[0], T2Bins2[-1] )
  353. ax1.set_ylim( ifaces[-1], ifaces[0] )
  354. ax2 = ax1.twiny()
  355. ax2.plot( np.sum(DINV, axis=1), (ifaces[1:]+ifaces[0:-1])/2 , color='mediumblue' )
  356. ax2.set_xlabel(u"total water (m$^3$/m$^3$)", color='mediumblue')
  357. ax2.set_ylim( ifaces[-1], ifaces[0] )
  358. ax2.xaxis.set_major_locator( MaxNLocator(nbins = 3) )
  359. ax2.get_xaxis().set_major_formatter(FormatStrFormatter('%0.2f'))
  360. pdf.savefig(facecolor=[0,0,0,0])
  361. plt.close(fig)
  362. np.save("MRM", MRM)
  363. centres = (ifaces[0:-1]+ifaces[1:])/2
  364. X,Y = np.meshgrid(ifaces,ifaces)
  365. fig = plt.figure( figsize=(pc2in(20.0),pc2in(22.)) )
  366. ax1 = fig.add_axes( [.2,.15,.6,.7] )
  367. ax1.pcolor(X,Y,MRM, cmap = cmocean.cm.ice)
  368. ax1.set_ylim(ifaces[-1], ifaces[0])
  369. maxDepth = np.argmax(MRM, axis=0)
  370. plt.plot(centres[maxDepth], centres, color='white')
  371. # Determine DOI
  372. DOIMetric = centres[maxDepth]/centres #> 0.9
  373. DOI = ifaces[ np.where(DOIMetric < 0.9 ) ][0]
  374. plt.axhline(y=DOI, color='white', linestyle='-.')
  375. ax1.set_ylim( ifaces[-1], ifaces[0] )
  376. ax1.set_xlim( ifaces[0], ifaces[-1] )
  377. ax1.set_xlabel(u"depth (m)")
  378. ax1.set_ylabel(u"depth (m)")
  379. plt.savefig("resolutionmatrix.pdf")
  380. pdf.close()
  381. INV = np.reshape(inv, (len(ifaces)-1,cont["T2Bins"]["number"]) )
  382. ############## LINEAR RESULT ##########################
  383. Y,X = meshgrid( ifaces, T2Bins2 )
  384. fig = plt.figure( figsize=(pc2in(20.0),pc2in(22.)) )
  385. ax1 = fig.add_axes( [.2,.15,.6,.7] )
  386. im = ax1.pcolor(X, Y, INV.T, cmap=cmocean.cm.tempo) #cmap='viridis')
  387. im.set_edgecolor('face')
  388. ax1.set_xlim( T2Bins[0], T2Bins[-1] )
  389. ax1.set_ylim( ifaces[-1], ifaces[0] )
  390. cb = plt.colorbar(im, label = u"PWC (m$^3$/m$^3$)") #, format='%1.1f')
  391. cb.locator = MaxNLocator( nbins = 4)
  392. cb.ax.yaxis.set_offset_position('left')
  393. cb.update_ticks()
  394. ax1.set_xlabel(u"$T_2^*$ (ms)")
  395. ax1.set_ylabel(u"depth (m)")
  396. ax1.get_xaxis().set_major_formatter(FormatStrFormatter('%1.0f'))
  397. ax1.get_yaxis().set_major_formatter(FormatStrFormatter('%1.0f'))
  398. ax1.xaxis.set_major_locator( MaxNLocator(nbins = 4) )
  399. ax2 = ax1.twiny()
  400. ax2.plot( np.sum(INV, axis=1), (ifaces[1:]+ifaces[0:-1])/2 , color='mediumblue' )
  401. ax2.set_xlabel(u"NMR total water (m$^3$/m$^3$)", color='mediumblue')
  402. ax2.set_ylim( ifaces[-1], ifaces[0] )
  403. ax2.xaxis.set_major_locator( MaxNLocator(nbins = 3) )
  404. ax2.get_xaxis().set_major_formatter(FormatStrFormatter('%0.2f'))
  405. #ax2.axhline( y=ifaces[SNRidx], xmin=0, xmax=1, color='black', linestyle='dashed' )
  406. if CalcDOI:
  407. ax2.axhline( y=DOI, xmin=0, xmax=1, color='black', linestyle='dashed' )
  408. ax2.tick_params(axis='x', colors='mediumblue')
  409. plt.setp(ax2.get_xticklabels(), color="mediumblue")
  410. plt.savefig("akvoInversion.pdf")
  411. #############
  412. # water plot#
  413. fig2 = plt.figure( figsize=(pc2in(20.0),pc2in(22.)) )
  414. ax = fig2.add_axes( [.2,.15,.6,.7] )
  415. # Bound water cutoff
  416. Bidx = T2Bins<33.0
  417. twater = np.sum(INV, axis=1)
  418. bwater = np.sum(INV[:,Bidx], axis=1)
  419. ax.plot( twater, (ifaces[0:-1]+ifaces[1::])/2, label="NMR total water", color='blue' )
  420. ax.plot( bwater, (ifaces[0:-1]+ifaces[1::])/2, label="NMR bound water", color='green' )
  421. ax.fill_betweenx((ifaces[0:-1]+ifaces[1::])/2 , twater, bwater, where=twater >= bwater, facecolor='blue', alpha=.5)
  422. ax.fill_betweenx((ifaces[0:-1]+ifaces[1::])/2 , bwater, 0, where=bwater >= 0, facecolor='green', alpha=.5)
  423. ax.set_xlabel(r"$\theta_N$ (m$^3$/m$^3$)")
  424. ax.set_ylabel(r"depth (m)")
  425. ax.set_ylim( ifaces[-1], ifaces[0] )
  426. ax.set_xlim( 0, ax.get_xlim()[1] )
  427. #ax.axhline( y=ifaces[SNRidx], xmin=0, xmax=1, color='black', linestyle='dashed' )
  428. if CalcDOI:
  429. ax.axhline( y=DOI, xmin=0, xmax=1, color='black', linestyle='dashed' )
  430. # Hide the right and top spines
  431. ax.spines['right'].set_visible(False)
  432. ax.spines['top'].set_visible(False)
  433. # Only show ticks on the left and bottom spines
  434. ax.yaxis.set_ticks_position('left')
  435. ax.xaxis.set_ticks_position('bottom')
  436. plt.savefig("akvoInversionWC.pdf")
  437. plt.legend()
  438. fr = pd.DataFrame( INV, columns=T2Bins ) #[0:-1] )
  439. fr.insert(0, "layer top", ifaces[0:-1] )
  440. fr.insert(1, "layer bottom", ifaces[1::] )
  441. fr.insert(2, "NMR total water", np.sum(INV, axis=1) )
  442. fr.insert(3, "NMR bound water", bwater )
  443. fr.insert(4, "Layer SNR", SNR )
  444. if CalcDOI:
  445. fr.insert(5, "Resolution", DOIMetric )
  446. fr.to_csv("akvoInversion.csv", mode='w+')
  447. ############## NONLINEAR RESULT ##########################
  448. if nonLinearRefinement:
  449. Y,X = meshgrid( ifaces, T2Bins2 )
  450. fig = plt.figure( figsize=(pc2in(20.0),pc2in(22.)) )
  451. ax1 = fig.add_axes( [.2,.15,.6,.7] )
  452. im = ax1.pcolor(X, Y, INVc.T, cmap=cmocean.cm.tempo) #cmap='viridis')
  453. #im = ax1.pcolor(X[0:SNRidx,:], Y[0:SNRidx,:], INV.T[0:SNRidx,:], cmap=cmocean.cm.tempo) #cmap='viridis')
  454. #im = ax1.pcolor(X[SNRidx::,:], Y[SNRidx::,:], INV.T[SNRidx::,:], cmap=cmocean.cm.tempo, alpha=.5) #cmap='viridis')
  455. #im = ax1.pcolormesh(X, Y, INV.T, alpha=alphas) #, cmap=cmocean.cm.tempo) #cmap='viridis')
  456. #im = ax1.pcolormesh(X, Y, INV.T, alpha=alphas) #, cmap=cmocean.cm.tempo) #cmap='viridis')
  457. #ax1.axhline( y=ifaces[SNRidx], xmin=T2Bins[0], xmax=T2Bins[-1], color='black' )
  458. im.set_edgecolor('face')
  459. ax1.set_xlim( T2Bins[0], T2Bins[-1] )
  460. ax1.set_ylim( ifaces[-1], ifaces[0] )
  461. cb = plt.colorbar(im, label = u"PWC (m$^3$/m$^3$)") #, format='%1.1f')
  462. cb.locator = MaxNLocator( nbins = 4)
  463. cb.ax.yaxis.set_offset_position('left')
  464. cb.update_ticks()
  465. ax1.set_xlabel(u"$T_2^*$ (ms)")
  466. ax1.set_ylabel(u"depth (m)")
  467. ax1.get_xaxis().set_major_formatter(FormatStrFormatter('%1.0f'))
  468. ax1.get_yaxis().set_major_formatter(FormatStrFormatter('%1.0f'))
  469. ax1.xaxis.set_major_locator( MaxNLocator(nbins = 4) )
  470. #ax1.xaxis.set_label_position('top')
  471. ax2 = ax1.twiny()
  472. ax2.plot( np.sum(INVc, axis=1), (ifaces[1:]+ifaces[0:-1])/2 , color='mediumblue' )
  473. ax2.set_xlabel(u"NMR total water (m$^3$/m$^3$)", color='mediumblue')
  474. ax2.set_ylim( ifaces[-1], ifaces[0] )
  475. ax2.xaxis.set_major_locator( MaxNLocator(nbins = 3) )
  476. ax2.get_xaxis().set_major_formatter(FormatStrFormatter('%0.2f'))
  477. #ax2.axhline( y=ifaces[SNRidx], xmin=0, xmax=1, color='black', linestyle='dashed' )
  478. if CalcDOI:
  479. ax2.axhline( y=DOI, xmin=0, xmax=1, color='black', linestyle='dashed' )
  480. #ax2.xaxis.set_label_position('bottom')
  481. #fig.suptitle("Non linear inversion")
  482. ax2.tick_params(axis='x', colors='mediumblue')
  483. plt.setp(ax2.get_xticklabels(), color="mediumblue")
  484. plt.savefig("akvoInversionNL.pdf")
  485. #############
  486. # water plot#
  487. fig2 = plt.figure( figsize=(pc2in(20.0),pc2in(22.)) )
  488. ax = fig2.add_axes( [.2,.15,.6,.7] )
  489. # Bound water cutoff
  490. Bidx = T2Bins<33.0
  491. twater = np.sum(INVc, axis=1)
  492. bwater = np.sum(INVc[:,Bidx], axis=1)
  493. ax.plot( twater, (ifaces[0:-1]+ifaces[1::])/2, label="NMR total water", color='blue' )
  494. ax.plot( bwater, (ifaces[0:-1]+ifaces[1::])/2, label="NMR bound water", color='green' )
  495. ax.fill_betweenx((ifaces[0:-1]+ifaces[1::])/2 , twater, bwater, where=twater >= bwater, facecolor='blue', alpha=.5)
  496. ax.fill_betweenx((ifaces[0:-1]+ifaces[1::])/2 , bwater, 0, where=bwater >= 0, facecolor='green', alpha=.5)
  497. ax.set_xlabel(r"$\theta_N$ (m$^3$/m$^3$)")
  498. ax.set_ylabel(r"depth (m)")
  499. ax.set_ylim( ifaces[-1], ifaces[0] )
  500. ax.set_xlim( 0, ax.get_xlim()[1] )
  501. #ax.axhline( y=ifaces[SNRidx], xmin=0, xmax=1, color='black', linestyle='dashed' )
  502. if CalcDOI:
  503. ax.axhline( y=DOI, xmin=0, xmax=1, color='black', linestyle='dashed' )
  504. # Hide the right and top spines
  505. ax.spines['right'].set_visible(False)
  506. ax.spines['top'].set_visible(False)
  507. # Only show ticks on the left and bottom spines
  508. ax.yaxis.set_ticks_position('left')
  509. ax.xaxis.set_ticks_position('bottom')
  510. plt.savefig("akvoNLInversionWC.pdf")
  511. plt.legend()
  512. # Report results into a text file
  513. fr = pd.DataFrame( INVc, columns=T2Bins ) #[0:-1] )
  514. fr.insert(0, "layer top", ifaces[0:-1] )
  515. fr.insert(1, "layer bottom", ifaces[1::] )
  516. fr.insert(2, "NMR total water", np.sum(INVc, axis=1) )
  517. fr.insert(3, "NMR bound water", bwater )
  518. fr.insert(4, "Layer SNR", SNR )
  519. if CalcDOI:
  520. fr.insert(5, "Resolution", DOIMetric )
  521. fr.to_csv("akvoNLInversion.csv", mode='w+')
  522. #fr.to_excel("akvoInversion.xlsx")
  523. plt.show()
  524. if __name__ == "__main__":
  525. main()