Browse Source

Added non-linear refinement to inversion, and improved inversion output

tags/1.6.1
Trevor Irons 3 years ago
parent
commit
b6316a1ebe
5 changed files with 171 additions and 23 deletions
  1. 2
    1
      akvo/tressel/calcAkvoKernel.py
  2. 103
    12
      akvo/tressel/invertTA.py
  3. 17
    6
      akvo/tressel/logbarrier.py
  4. 45
    0
      akvo/tressel/nonlinearinv.py
  5. 4
    4
      setup.py

+ 2
- 1
akvo/tressel/calcAkvoKernel.py View File

@@ -64,6 +64,7 @@ def main():
64 64
 
65 65
     if len(sys.argv) < 2:
66 66
         print ("usage  python calcAkvoKernel.py   AkvoDataset.yaml  Coil1.yaml kparams.yaml  SaveString.yaml " )
67
+        print ("usage  akvoKO   AkvoDataset.yaml   kparams.yaml  SaveString.yaml " )
67 68
         exit()
68 69
 
69 70
     AKVO = loadAkvoData(sys.argv[1])
@@ -100,7 +101,7 @@ def main():
100 101
             Coil1.SetCurrent(1.)
101 102
             Kern.PushCoil( rx.split('.yml')[0], Coil1 )
102 103
         else:
103
-            print("reuse")
104
+            print("reuse tx coil")
104 105
         RX.append( rx.split('.yml')[0] )
105 106
     
106 107
 

+ 103
- 12
akvo/tressel/invertTA.py View File

@@ -17,7 +17,8 @@ from matplotlib.ticker import FormatStrFormatter
17 17
 from matplotlib.colors import Normalize
18 18
 
19 19
 import cmocean
20
-from akvo.tressel.lemma_yaml import * 
20
+from akvo.tressel.lemma_yaml import *
21
+from akvo.tressel import nonlinearinv as nl 
21 22
 
22 23
 import pandas as pd
23 24
 
@@ -28,7 +29,7 @@ def buildKQT(K0,tg,T2Bins):
28 29
     nlay, nq = np.shape(K0)
29 30
     nt2 = len(T2Bins)
30 31
     nt = len(tg)
31
-    KQT = np.zeros( ( nq*nt,nt2*nlay) )
32
+    KQT = np.zeros( ( nq*nt,nt2*nlay), dtype=np.complex128 )
32 33
     for iq in range(nq):
33 34
         for it in range(nt):
34 35
             for ilay in range(nlay):
@@ -87,7 +88,8 @@ def loadK0(fname):
87 88
         K0 = yaml.load(f, Loader=yaml.Loader)
88 89
     K = catLayers(K0.K0)
89 90
     ifaces = np.array(K0.Interfaces.data)
90
-    return ifaces, np.abs(K)
91
+    return ifaces, K
92
+    #return ifaces, np.abs(K)
91 93
 
92 94
 
93 95
 
@@ -131,12 +133,15 @@ def main():
131 133
     for ik in range(1, len(K0)):
132 134
         K0[0] = np.concatenate( (K0[0].T, K0[ik].T) ).T
133 135
     K0 = K0[0]
134
-    #plt.matshow(K0)
136
+
137
+    #plt.matshow(np.real(K0))
138
+    #plt.show()
139
+    #exit()
135 140
 
136 141
     ###################    
137 142
     # VERY Simple DOI #
138
-    maxq = np.argmax(K0, axis=1)
139
-    maxK = .1 *  K0[ np.arange(0,len(ifaces)-1), maxq ] # 10% water is arbitrary  
143
+    maxq = np.argmax(np.abs(K0), axis=1)
144
+    maxK = .1 *  np.abs(K0)[ np.arange(0,len(ifaces)-1), maxq ] # 10% water is arbitrary  
140 145
     SNR = maxK / (VS[0][0])
141 146
 
142 147
     #SNR[SNR>1] = 1
@@ -163,18 +168,60 @@ def main():
163 168
     # Build full kernel
164 169
     ############################################### 
165 170
     T2Bins = np.logspace( np.log10(cont["T2Bins"]["low"]), np.log10(cont["T2Bins"]["high"]), cont["T2Bins"]["number"], endpoint=True, base=10)  
166
-    KQT = buildKQT(K0,tg,T2Bins)
167
- 
171
+    KQT = np.real(buildKQT(np.abs(K0),tg,T2Bins))
172
+
173
+    # model resolution matrix 
174
+    np.linalg.svd(KQT)
175
+
176
+    exit()
177
+
168 178
     ###############################################
169
-    # Invert
179
+    # Linear Inversion 
170 180
     ############################################### 
171 181
     print("Calling inversion", flush=True)
172
-    inv, ibreak, errn, phim, phid, mkappa = logBarrier(KQT, np.ravel(V), T2Bins, "lcurve", MAXITER=150, sigma=np.ravel(VS), alpha=1e6, smooth="Smallest" ) 
182
+    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" ) 
173 183
 
184
+    ###############################################
185
+    # Non-linear refinement! 
186
+    ###############################################   
187
+ 
188
+    KQTc = buildKQT(K0, tg, T2Bins)
189
+    prec = np.abs(np.dot(KQTc, inv))
190
+    phidc = np.linalg.norm(np.dot(Wd,prec-np.ravel(V)))**2
191
+    #PREc = np.reshape( prec, np.shape(V)  )
192
+    print("PHID linear=", errn, "PHID complex=", phidc/len(np.ravel(V)))
193
+    
194
+    res = nl.nonlinearinversion(inv, Wd, KQTc, np.ravel(V), Wm, alphastar )   
195
+    if res.success == True:    
196
+        INVc = np.reshape(res.x, (len(ifaces)-1,cont["T2Bins"]["number"]) )
197
+        prec = np.abs(np.dot(KQTc, res.x))
198
+        phidc = np.linalg.norm(np.dot(Wd,prec-np.ravel(V)))**2
199
+        #PREc = np.reshape( prec, np.shape(V)  )
200
+        print("PHID linear=", errn, "PHID nonlinear=", phidc/len(np.ravel(V)))
201
+   
202
+    # Perform second fit around results of first  
203
+    res = nl.nonlinearinversion(res.x, Wd, KQTc, np.ravel(V), Wm, alphastar )   
204
+    if res.success == True:    
205
+        INVc = np.reshape(res.x, (len(ifaces)-1,cont["T2Bins"]["number"]) )
206
+        prec = np.abs(np.dot(KQTc, res.x))
207
+        phidc = np.linalg.norm(np.dot(Wd,prec-np.ravel(V)))**2
208
+        #PREc = np.reshape( prec, np.shape(V)  )
209
+        print("PHID linear=", errn, "PHID nonlinear=", phidc/len(np.ravel(V)))
210
+
211
+    #plt.matshow(INVc)
212
+    #KQTc = buildKQT(K0,tg,T2Bins)
213
+
214
+    #plt.matshow(PREc, cmap='Blues')
215
+    #plt.gca().set_title("complex predicted")
216
+    #plt.colorbar()
174 217
 
175 218
     ###############################################
176 219
     # Appraise
177 220
     ###############################################
221
+
222
+
223
+
224
+
178 225
  
179 226
     pre = np.dot(KQT,inv) 
180 227
     PRE = np.reshape( pre, np.shape(V)  )
@@ -192,8 +239,8 @@ def main():
192 239
     plt.gca().set_title("observed")
193 240
     plt.colorbar()
194 241
 
242
+
195 243
     T2Bins = np.append( T2Bins, T2Bins[-1] + (T2Bins[-1]-T2Bins[-2]) )
196
-    
197 244
     INV = np.reshape(inv, (len(ifaces)-1,cont["T2Bins"]["number"]) )
198 245
 
199 246
     #alphas = np.tile(SNR, (len(T2Bins)-1,1))
@@ -205,6 +252,8 @@ def main():
205 252
 
206 253
     #greys = np.full((*(INV.T).shape, 3), 70, dtype=np.uint8)
207 254
 
255
+    ##############  LINEAR RESULT   ##########################
256
+
208 257
     Y,X = meshgrid( ifaces, T2Bins )
209 258
     fig = plt.figure( figsize=(pc2in(20.0),pc2in(22.)) )
210 259
     ax1 = fig.add_axes( [.2,.15,.6,.7] )
@@ -242,6 +291,48 @@ def main():
242 291
 
243 292
     plt.savefig("akvoInversion.pdf")
244 293
 
294
+
295
+    ##############  NONLINEAR RESULT   ##########################
296
+
297
+    Y,X = meshgrid( ifaces, T2Bins )
298
+    fig = plt.figure( figsize=(pc2in(20.0),pc2in(22.)) )
299
+    ax1 = fig.add_axes( [.2,.15,.6,.7] )
300
+    im = ax1.pcolor(X, Y, INVc.T, cmap=cmocean.cm.tempo) #cmap='viridis')
301
+    #im = ax1.pcolor(X[0:SNRidx,:], Y[0:SNRidx,:], INV.T[0:SNRidx,:], cmap=cmocean.cm.tempo) #cmap='viridis')
302
+    #im = ax1.pcolor(X[SNRidx::,:], Y[SNRidx::,:], INV.T[SNRidx::,:], cmap=cmocean.cm.tempo, alpha=.5) #cmap='viridis')
303
+    #im = ax1.pcolormesh(X, Y, INV.T, alpha=alphas) #, cmap=cmocean.cm.tempo) #cmap='viridis')
304
+    #im = ax1.pcolormesh(X, Y, INV.T, alpha=alphas) #, cmap=cmocean.cm.tempo) #cmap='viridis')
305
+    #ax1.axhline( y=ifaces[SNRidx], xmin=T2Bins[0], xmax=T2Bins[-1], color='black'  )
306
+    im.set_edgecolor('face')
307
+    ax1.set_xlim( T2Bins[0], T2Bins[-1] )
308
+    ax1.set_ylim( ifaces[-1], ifaces[0] )
309
+    cb = plt.colorbar(im, label = u"PWC (m$^3$/m$^3$)") #, format='%1.1f')
310
+    cb.locator = MaxNLocator( nbins = 4)
311
+    cb.ax.yaxis.set_offset_position('left')                         
312
+    cb.update_ticks()
313
+ 
314
+    ax1.set_xlabel(u"$T_2^*$ (ms)")
315
+    ax1.set_ylabel(u"depth (m)")
316
+    
317
+    ax1.get_xaxis().set_major_formatter(FormatStrFormatter('%1.0f'))
318
+    ax1.get_yaxis().set_major_formatter(FormatStrFormatter('%1.0f'))
319
+    ax1.xaxis.set_major_locator( MaxNLocator(nbins = 4) )   
320
+
321
+    #ax1.xaxis.set_label_position('top') 
322
+
323
+    ax2 = ax1.twiny()
324
+    ax2.plot( np.sum(INVc, axis=1), (ifaces[1:]+ifaces[0:-1])/2 ,  color='red' )
325
+    ax2.set_xlabel(u"total water (m$^3$/m$^3$)")
326
+    ax2.set_ylim( ifaces[-1], ifaces[0] )
327
+    ax2.xaxis.set_major_locator( MaxNLocator(nbins = 3) )   
328
+    ax2.get_xaxis().set_major_formatter(FormatStrFormatter('%0.2f'))
329
+    #ax2.axhline( y=ifaces[SNRidx], xmin=0, xmax=1, color='black', linestyle='dashed'  )
330
+    #ax2.xaxis.set_label_position('bottom') 
331
+    fig.suptitle("Non linear inversion")
332
+    plt.savefig("akvoInversionNL.pdf")
333
+
334
+
335
+
245 336
     #############
246 337
     # water plot#
247 338
 
@@ -265,7 +356,7 @@ def main():
265 356
     ax.set_ylim( ifaces[-1], ifaces[0] )
266 357
     ax.set_xlim( 0, ax.get_xlim()[1] )
267 358
     
268
-    ax.axhline( y=ifaces[SNRidx], xmin=0, xmax=1, color='black', linestyle='dashed'  )
359
+    #ax.axhline( y=ifaces[SNRidx], xmin=0, xmax=1, color='black', linestyle='dashed'  )
269 360
     
270 361
     plt.savefig("akvoInversionWC.pdf")
271 362
     plt.legend()  

+ 17
- 6
akvo/tressel/logbarrier.py View File

@@ -128,6 +128,8 @@ def logBarrier(A, b, T2Bins, lambdastar, x_0=0, xr=0, alpha=10, mu1=10, mu2=10,
128 128
     ALPHA = []
129 129
     ALPHA.append(alpha)
130 130
     #ALPHA = np.linspace( alpha, 1, MAXITER  )
131
+    print ("{:^10} {:^15} {:^15} {:^15} {:^15} {:^10} {:^10}".format("iteration",  "lambda", "phi_d", "phi_m","phi","kappa","kappa dist."), flush=True) 
132
+    print ("{:^10} {:>15} {:<15} {:<15} {:<15} {:<10} {:<10}".format("----------", "---------------", "---------------","---------------","---------------","----------","----------"), flush=True) 
131 133
     for i in range(MAXITER):
132 134
         #alpha = ALPHA[i]
133 135
 
@@ -202,6 +204,8 @@ def logBarrier(A, b, T2Bins, lambdastar, x_0=0, xr=0, alpha=10, mu1=10, mu2=10,
202 204
         PHID.append(phid)      
203 205
         MOD.append(np.copy(x))  
204 206
 
207
+        tphi = phid + alpha*phim
208
+
205 209
         # determine alpha
206 210
         scale = 1.5*(len(b)/phid)
207 211
         #alpha *= np.sqrt(scale)
@@ -211,8 +215,11 @@ def logBarrier(A, b, T2Bins, lambdastar, x_0=0, xr=0, alpha=10, mu1=10, mu2=10,
211 215
         ALPHA.append(alpha)
212 216
         #alpha = ALPHA[i+1]
213 217
             
214
-        print("inversion progress", i, alpha, np.sqrt(phid/len(b)), phim, flush=True)       
218
+        #print("inversion progress", i, alpha, np.sqrt(phid/len(b)), phim, flush=True)      
219
+        #print ("{:<8} {:<15} {:<10} {:<10}".format(i, alpha, np.sqrt(phid/len(b)), phim), flush=True) 
215 220
         
221
+        if i < 4:        
222
+            print ("{:^10} {:>15.4f} {:>15.4f} {:>15.4f} {:>15.4f}".format(i, alpha, np.sqrt(phid/len(b)), phim, tphi ), flush=True) 
216 223
 
217 224
 #         if np.sqrt(phid/len(b)) < 0.97: 
218 225
 #             ibreak = -1
@@ -234,23 +241,25 @@ def logBarrier(A, b, T2Bins, lambdastar, x_0=0, xr=0, alpha=10, mu1=10, mu2=10,
234 241
             if i > 4: 
235 242
                 kappa = curvaturefd(np.log(np.array(PHIM)), np.log(np.array(PHID)), ALPHA[0:i+1])#ALPHA[0:-1])
236 243
                 #kappa = curvatureg(np.log(np.array(PHIM)), np.log(np.array(PHID)))
237
-                print("max kappa", np.argmax(kappa), "distance from", i-np.argmax(kappa)) 
244
+                #print("max kappa", np.argmax(kappa), "distance from", i-np.argmax(kappa)) 
245
+                print ("{:^10} {:>15.4f} {:>15.4f} {:>15.4f} {:>15.4f} {:^10} {:^10}".format(i, alpha, np.sqrt(phid/len(b)), phim, tphi, np.argmax(kappa), i-np.argmax(kappa)), flush=True) 
238 246
             if i > 4 and (i-np.argmax(kappa)) > 4: # ((np.sqrt(phid_old/len(b))-np.sqrt(phid/len(b))) < 1e-4) : 
239 247
             #if np.sqrt(phid/len(b)) < 3.0 and ((np.sqrt(phid_old/len(b))-np.sqrt(phid/len(b))) < 1e-3): 
240 248
                 ibreak = 1
241 249
                 MOD = np.array(MOD)
242
-                print ("###########################") #slow convergence", alpha, "phid_old", np.sqrt(phid_old/len(b)), "phid", np.sqrt(phid/len(b)), ibreak)
250
+                print ("################################") #slow convergence", alpha, "phid_old", np.sqrt(phid_old/len(b)), "phid", np.sqrt(phid/len(b)), ibreak)
243 251
                 print ("Using L-curve criteria") 
244 252
                 #kappa = curvaturefd(np.log(np.array(PHIM)), np.log(np.array(PHID)), ALPHA[0:-1])
245 253
                 #kappa2 = curvatureg(np.log(np.array(PHIM)), np.log(np.array(PHID)))
246 254
                 #kappa = curvature( np.array(PHIM), np.array(PHID))
247 255
                 x = MOD[ np.argmax(kappa) ]
256
+                alphastar = ALPHA[ np.argmax(kappa) ]
248 257
                 b_pre = np.dot(A, x)
249 258
                 phid = np.linalg.norm( np.dot(Wd, (b-b_pre)))**2
250 259
                 phim = np.linalg.norm( np.dot(Phim_base, (x-xr)) )**2
251 260
                 mu1 = ((phid + alpha*phim) / phib) 
252
-                print ("L-curve selected", alpha, "phid_old", np.sqrt(phid_old/len(b)), "phid", np.sqrt(phid/len(b)), ibreak)
253
-                print ("###########################")
261
+                print ("L-curve selected: iteration=", np.argmax(kappa)) #, " lambda*=", alpha, "phid_old=", np.sqrt(phid_old/len(b)), "phid=", np.sqrt(phid/len(b)), ibreak)
262
+                print ("################################")
254 263
                 if np.sqrt(phid/len(b)) <= 1:
255 264
                     ibreak=0
256 265
 
@@ -291,8 +300,10 @@ def logBarrier(A, b, T2Bins, lambdastar, x_0=0, xr=0, alpha=10, mu1=10, mu2=10,
291 300
         mu1 = ((phid + alpha*phim) / phib) 
292 301
 
293 302
     if lambdastar == "lcurve":
294
-        return x, ibreak, np.sqrt(phid/len(b)), PHIM, PHID/len(b), np.argmax(kappa)
303
+        #print("Returning L curve result")
304
+        return x, ibreak, np.sqrt(phid/len(b)), PHIM, PHID/len(b), np.argmax(kappa), Wd, Phim_base, alphastar
295 305
     else:
306
+        print("")
296 307
         return x, ibreak, np.sqrt(phid/len(b))
297 308
 
298 309
 

+ 45
- 0
akvo/tressel/nonlinearinv.py View File

@@ -0,0 +1,45 @@
1
+import numpy as np
2
+import scipy.optimize as so 
3
+
4
+def phid(Wd, K, m, d_obs):
5
+    """
6
+        Wd = data weighting matrix 
7
+        K = complex valued forward model kernel 
8
+        m = model 
9
+        d_obs = observed data 
10
+    """
11
+    #print("phid=", np.linalg.norm(np.dot(Wd, np.abs(np.dot(K,m)) - d_obs))**2 / len(d_obs) )
12
+    return np.linalg.norm(np.dot(Wd, np.abs(np.dot(K,m)) - d_obs))**2
13
+
14
+def phim(Wm, m):
15
+    """
16
+        Wm = model weighting matrix 
17
+        x = model 
18
+    """
19
+    return np.linalg.norm(np.dot(Wm, m))**2
20
+
21
+def PHI(m, Wd, K, d_obs, Wm, alphastar):
22
+    """
23
+        Global objective function 
24
+        x = model to be fit 
25
+        Wd = data weighting matrix 
26
+        K = complex forward modelling kernel 
27
+        d_obs = observed data (modulus)
28
+        Wm = model weighting matrix 
29
+        alphastar = regularisation to use
30
+    """
31
+    return phid(Wd, K, m, d_obs) + alphastar*phim(Wm, m)
32
+
33
+# main 
34
+def nonlinearinversion( x0, Wd, K, d_obs, Wm, alphastar):
35
+    print("Performing non-linear inversion")
36
+    args = (Wd, K, d_obs, Wm, alphastar)
37
+    #return so.minimize(PHI, np.zeros(len(x0)), args, 'Nelder-Mead')
38
+    bounds = np.zeros((len(x0),2))
39
+    bounds[:,0] = x0*0.75
40
+    bounds[:,1] = x0*1.25
41
+    return so.minimize(PHI, x0, args, 'L-BFGS-B', bounds=bounds)     # Works well 
42
+    #return so.minimize(PHI, x0, args, 'Powell', bounds=bounds)       # Slow but works 
43
+    #return so.minimize(PHI, x0, args, 'trust-constr', bounds=bounds) # very Slow 
44
+    #return so.minimize(PHI, x0, args, 'TNC', bounds=bounds)          # slow 
45
+    #return so.minimize(PHI, x0, args, 'SLSQP', bounds=bounds)        # slow 

+ 4
- 4
setup.py View File

@@ -20,9 +20,9 @@ class custom_build_py(build_py):
20 20
 with open("README.md", "r") as fh:
21 21
     long_description = fh.read()
22 22
 
23
-setup(name='Akvo',
24
-      version='1.5.2',
25
-      python_requires='>3.7.0', # due to pyLemma
23
+setup(name='Akvo',     
24
+      version='1.5.2', 
25
+      python_requires='>3.7.0', # due to pyLemma 
26 26
       description='Surface nuclear magnetic resonance workbench',
27 27
       long_description=long_description,
28 28
       long_description_content_type='text/markdown',
@@ -56,7 +56,7 @@ setup(name='Akvo',
56 56
           'pyLemma >= 0.4.0'
57 57
       ],
58 58
       packages=['akvo', 'akvo.tressel', 'akvo.gui'],
59
-      license=['GPL 4.0'],
59
+      license='GPL 4.0',
60 60
       entry_points = {
61 61
               'console_scripts': [
62 62
                   'akvo = akvo.gui.akvoGUI:main',                  

Loading…
Cancel
Save