Surface NMR processing and inversion GUI
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.

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165
  1. import numpy as np
  2. from scipy.optimize import least_squares
  3. from scipy.optimize import minimize
  4. from scipy.linalg import lstsq as sclstsq
  5. import scipy.linalg as lin
  6. def harmonicEuler ( sN, fs, t, f0, k1, kN, ks ):
  7. """
  8. Performs inverse calculation of harmonics contaminating a signal.
  9. Args:
  10. sN = signal containing noise
  11. fs = sampling frequency
  12. t = time samples
  13. f0 = base frequency of the sinusoidal noise
  14. nK = number of harmonics to calculate
  15. """
  16. KK = np.arange(k1, kN+1, 1/ks )
  17. nK = len(KK)
  18. A = np.exp(1j* np.tile(KK,(len(t), 1)) * 2*np.pi* (f0/fs) * np.tile(np.arange(1, len(t)+1, 1),(nK,1)).T)
  19. v = np.linalg.lstsq(A, sN, rcond=None)
  20. alpha = np.real(v[0])
  21. beta = np.imag(v[0])
  22. amp = np.abs(v[0])
  23. phase = np.angle(v[0])
  24. h = np.zeros(len(t))
  25. for ik, k in enumerate(KK):
  26. h += 2*amp[ik] * np.cos( 2.*np.pi*(k) * (f0/fs) * np.arange(1, len(t)+1, 1 ) + phase[ik] )
  27. return sN-h
  28. def harmonicNorm (f0, sN, fs, t, k1, kN, ks):
  29. #return np.linalg.norm( harmonicEuler(sN, fs, t, f0, k1, kN, ks))
  30. ii = sN < (3.* np.std(sN))
  31. return np.linalg.norm( harmonicEuler(sN, fs, t, f0, k1, kN, ks)[ii] )
  32. def minHarmonic(sN, fs, t, f0, k1, kN, ks):
  33. # CG, BFGS, Newton-CG, L-BFGS-B, TNC, SLSQP, dogleg, trust-ncg, trust-krylov, trust-exact and trust-constr
  34. res = minimize(harmonicNorm, np.array((f0)), args=(sN, fs, t, k1, kN, ks), jac='2-point', method='BFGS') # hess=None, bounds=None )
  35. #print(res)
  36. #print ( "guess", guessf0( harmonicEuler(sN, fs, t, res.x[0], k1, kN, ks), fs ) )
  37. return harmonicEuler(sN, fs, t, res.x[0], k1, kN, ks)#[0]
  38. def harmonicEuler2 ( sN, fs, t, f0, f0k1, f0kN, f0ks, f1, f1k1, f1kN, f1ks ):
  39. """
  40. Performs inverse calculation of harmonics contaminating a signal.
  41. Args:
  42. sN = signal containing noise
  43. fs = sampling frequency
  44. t = time samples
  45. f0 = first base frequency of the sinusoidal noise
  46. f0k1 = First harmonic to calulate for f0
  47. f0kN = Last base harmonic to calulate for f0
  48. f0ks = subharmonics to calculate
  49. f1 = second base frequency of the sinusoidal noise
  50. f1k1 = First harmonic to calulate for f1
  51. f1kN = Last base harmonic to calulate for f1
  52. f1ks = subharmonics to calculate at f1 base frequency
  53. """
  54. KK0 = np.arange(f0k1, f0kN+1, 1/f0ks)
  55. nK0 = len(KK0)
  56. A0 = np.exp(1j* np.tile(KK0,(len(t), 1)) * 2*np.pi* (f0/fs) * np.tile( np.arange(1, len(t)+1, 1), (nK0,1)).T)
  57. KK1 = np.arange(f1k1, f1kN+1, 1/f1ks)
  58. nK1 = len(KK1)
  59. A1 = np.exp(1j* np.tile(KK1,(len(t), 1)) * 2*np.pi* (f1/fs) * np.tile( np.arange(1, len(t)+1, 1),(nK1,1)).T)
  60. A = np.concatenate((A0, A1), axis=1)
  61. v = np.linalg.lstsq(A, sN, rcond=None) # rcond=None) #, rcond=1e-8)
  62. amp0 = np.abs(v[0][0:nK0])
  63. phase0 = np.angle(v[0][0:nK0])
  64. amp1 = np.abs(v[0][nK0::])
  65. phase1 = np.angle(v[0][nK0::])
  66. h = np.zeros(len(t))
  67. for ik, k in enumerate(KK0):
  68. h += 2*amp0[ik] * np.cos( 2.*np.pi*(k) * (f0/fs) * np.arange(1, len(t)+1, 1 ) + phase0[ik] )
  69. for ik, k in enumerate(KK1):
  70. h += 2*amp1[ik] * np.cos( 2.*np.pi*(k) * (f0/fs) * np.arange(1, len(t)+1, 1 ) + phase1[ik] )
  71. return sN-h
  72. def harmonic2Norm (f0, sN, fs, t, f0k1, f0kN, f0ks, f1k1, f1kN, f1ks):
  73. #return np.linalg.norm(harmonicEuler2(f0[0], f0[1], sN, fs, nK, t))
  74. ii = sN < (3.* np.std(sN))
  75. return np.linalg.norm( harmonicEuler2(sN, fs, t, f0[0], f0k1, f0kN, f0ks, f0[1], f1k1, f1kN, f1ks)[ii] )
  76. def minHarmonic2(sN, fs, t, f0, f0k1, f0kN, f0ks, f1, f1k1, f1kN, f1ks):
  77. # CG, BFGS, Newton-CG, L-BFGS-B, TNC, SLSQP, dogleg, trust-ncg, trust-krylov, trust-exact and trust-constr
  78. res = minimize(harmonic2Norm, np.array((f0, f1)), args=(sN, fs, t, f0k1, f0kN, f0ks, f1k1,f1kN, f1ks), jac='2-point', method='BFGS') # hess=None, bounds=None )
  79. #print(res)
  80. #print ( "guess", guessf0(harmonicEuler2(sN, fs, t, res.x[0], f0k1, f0kN, f0ks, res.x[1], f1k1, f1kN, f1ks), fs) )
  81. return harmonicEuler2(sN, fs, t, res.x[0], f0k1, f0kN, f0ks, res.x[1], f1k1, f1kN, f1ks)#[0]
  82. def guessf0( sN, fs ):
  83. S = np.fft.fft(sN)
  84. w = np.fft.fftfreq( len(sN), 1/fs )
  85. imax = np.argmax( np.abs(S) )
  86. #plt.plot( w, np.abs(S) )
  87. #plt.show()
  88. #print(w)
  89. #print ( w[imax], w[imax+1] )
  90. return abs(w[imax])
  91. if __name__ == "__main__":
  92. import matplotlib.pyplot as plt
  93. f0 = 60 # Hz
  94. f1 = 60 # Hz
  95. delta = np.random.rand() - .5
  96. delta2 = np.random.rand() - .5
  97. print("delta", delta)
  98. print("delta2", delta2)
  99. fs = 10000 # GMR
  100. t = np.arange(0, 1, 1/fs)
  101. phi = 2.*np.pi*np.random.rand() - np.pi
  102. phi2 = 2.*np.pi*np.random.rand() - np.pi
  103. print("phi", phi, phi2)
  104. A = 1.0
  105. A2 = 0.0
  106. A3 = 1.0
  107. nK = 10
  108. T2 = .200
  109. sN = A *np.sin( ( 1*(delta +f0))*2*np.pi*t + phi ) + \
  110. A2*np.sin( ( 1*(delta2 +f1))*2*np.pi*t + phi2 ) + \
  111. np.random.normal(0,.1,len(t)) + \
  112. + A3*np.exp( -t/T2 )
  113. sNc = A *np.sin( (1*(delta +f0))*2*np.pi*t + phi ) + \
  114. A2*np.sin( (1*(delta2+f1))*2*np.pi*t + phi2 ) + \
  115. + A3*np.exp( -t/T2 )
  116. guessf0(sN, fs)
  117. # single freq
  118. #h = harmonicEuler( f0, sN, fs, nK, t)
  119. h = minHarmonic( f0, sN, fs, nK, t)
  120. # two freqs
  121. #h = minHarmonic2( f0+1e-2, f1-1e-2, sN, fs, nK, t)
  122. #h = harmonicEuler2( f0, f1, sN, fs, nK, t)
  123. plt.figure()
  124. plt.plot(t, sN, label="sN")
  125. #plt.plot(t, sN-h, label="sN-h")
  126. plt.plot(t, h, label='h')
  127. plt.title("harmonic")
  128. plt.legend()
  129. plt.figure()
  130. plt.plot(t, sN-sNc, label='true noise')
  131. plt.plot(t, h, label='harmonic removal')
  132. plt.plot(t, np.exp(-t/T2), label="nmr")
  133. plt.legend()
  134. plt.title("true noise")
  135. plt.show()