Surface NMR processing and inversion GUI
Vous ne pouvez pas sélectionner plus de 25 sujets Les noms de sujets doivent commencer par une lettre ou un nombre, peuvent contenir des tirets ('-') et peuvent comporter jusqu'à 35 caractères.

harmonic.py 6.6KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187
  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. #############
  36. # Reduced N #
  37. #############
  38. #print("min 10")
  39. #res = minimize(harmonicNorm, np.array((f0)), args=(sN, fs, t, k1, 10, ks), jac='2-point', method='BFGS') # hess=None, bounds=None )
  40. ##########
  41. # Bounds #
  42. ##########
  43. print("bounded search for f0", f0)
  44. bnds = ( (f0-0.125, f0+0.125), )
  45. #print("len f0", len( [f0,] ))
  46. #print("len bnds", len(bnds) )
  47. res = minimize(harmonicNorm, (f0,), args=(sN, fs, t, k1, kN, ks), jac='2-point', method='L-BFGS-B', bounds=bnds ) # hess=None, bounds=None )
  48. #print(res)
  49. #print ( "guess", guessf0( harmonicEuler(sN, fs, t, res.x[0], k1, kN, ks), fs ) )
  50. return harmonicEuler(sN, fs, t, res.x[0], k1, kN, ks), res.x[0]#[0]
  51. def harmonicEuler2 ( sN, fs, t, f0, f0k1, f0kN, f0ks, f1, f1k1, f1kN, f1ks ):
  52. """
  53. Performs inverse calculation of harmonics contaminating a signal.
  54. Args:
  55. sN = signal containing noise
  56. fs = sampling frequency
  57. t = time samples
  58. f0 = first base frequency of the sinusoidal noise
  59. f0k1 = First harmonic to calulate for f0
  60. f0kN = Last base harmonic to calulate for f0
  61. f0ks = subharmonics to calculate
  62. f1 = second base frequency of the sinusoidal noise
  63. f1k1 = First harmonic to calulate for f1
  64. f1kN = Last base harmonic to calulate for f1
  65. f1ks = subharmonics to calculate at f1 base frequency
  66. """
  67. KK0 = np.arange(f0k1, f0kN+1, 1/f0ks)
  68. nK0 = len(KK0)
  69. 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)
  70. KK1 = np.arange(f1k1, f1kN+1, 1/f1ks)
  71. nK1 = len(KK1)
  72. 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)
  73. A = np.concatenate((A0, A1), axis=1)
  74. v = np.linalg.lstsq(A, sN, rcond=None) # rcond=None) #, rcond=1e-8)
  75. amp0 = np.abs(v[0][0:nK0])
  76. phase0 = np.angle(v[0][0:nK0])
  77. amp1 = np.abs(v[0][nK0::])
  78. phase1 = np.angle(v[0][nK0::])
  79. h = np.zeros(len(t))
  80. for ik, k in enumerate(KK0):
  81. h += 2*amp0[ik] * np.cos( 2.*np.pi*(k) * (f0/fs) * np.arange(1, len(t)+1, 1 ) + phase0[ik] )
  82. for ik, k in enumerate(KK1):
  83. h += 2*amp1[ik] * np.cos( 2.*np.pi*(k) * (f0/fs) * np.arange(1, len(t)+1, 1 ) + phase1[ik] )
  84. return sN-h
  85. def harmonic2Norm (f0, sN, fs, t, f0k1, f0kN, f0ks, f1k1, f1kN, f1ks):
  86. #return np.linalg.norm(harmonicEuler2(f0[0], f0[1], sN, fs, nK, t))
  87. ii = sN < (3.* np.std(sN))
  88. return np.linalg.norm( harmonicEuler2(sN, fs, t, f0[0], f0k1, f0kN, f0ks, f0[1], f1k1, f1kN, f1ks)[ii] )
  89. def minHarmonic2(sN, fs, t, f0, f0k1, f0kN, f0ks, f1, f1k1, f1kN, f1ks):
  90. # CG, BFGS, Newton-CG, L-BFGS-B, TNC, SLSQP, dogleg, trust-ncg, trust-krylov, trust-exact and trust-constr
  91. #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 )
  92. print("search 10")
  93. res = minimize(harmonic2Norm, np.array((f0, f1)), args=(sN, fs, t, f0k1, 10, f0ks, f1k1, 10, f1ks), jac='2-point', method='BFGS') # hess=None, bounds=None )
  94. #print(res)
  95. #print ( "guess", guessf0(harmonicEuler2(sN, fs, t, res.x[0], f0k1, f0kN, f0ks, res.x[1], f1k1, f1kN, f1ks), fs) )
  96. return harmonicEuler2(sN, fs, t, res.x[0], f0k1, f0kN, f0ks, res.x[1], f1k1, f1kN, f1ks), res.x[0], res.x[1]#[0]
  97. def guessf0( sN, fs ):
  98. S = np.fft.fft(sN)
  99. w = np.fft.fftfreq( len(sN), 1/fs )
  100. imax = np.argmax( np.abs(S) )
  101. #plt.plot( w, np.abs(S) )
  102. #plt.show()
  103. #print(w)
  104. #print ( w[imax], w[imax+1] )
  105. return abs(w[imax])
  106. if __name__ == "__main__":
  107. import matplotlib.pyplot as plt
  108. f0 = 60 # Hz
  109. f1 = 60 # Hz
  110. delta = np.random.rand() - .5
  111. delta2 = np.random.rand() - .5
  112. print("delta", delta)
  113. print("delta2", delta2)
  114. fs = 10000 # GMR
  115. t = np.arange(0, 1, 1/fs)
  116. phi = 2.*np.pi*np.random.rand() - np.pi
  117. phi2 = 2.*np.pi*np.random.rand() - np.pi
  118. print("phi", phi, phi2)
  119. A = 1.0
  120. A2 = 0.0
  121. A3 = 1.0
  122. nK = 10
  123. T2 = .200
  124. sN = A *np.sin( ( 1*(delta +f0))*2*np.pi*t + phi ) + \
  125. A2*np.sin( ( 1*(delta2 +f1))*2*np.pi*t + phi2 ) + \
  126. np.random.normal(0,.1,len(t)) + \
  127. + A3*np.exp( -t/T2 )
  128. sNc = A *np.sin( (1*(delta +f0))*2*np.pi*t + phi ) + \
  129. A2*np.sin( (1*(delta2+f1))*2*np.pi*t + phi2 ) + \
  130. + A3*np.exp( -t/T2 )
  131. guessf0(sN, fs)
  132. # single freq
  133. #h = harmonicEuler( f0, sN, fs, nK, t)
  134. h = minHarmonic( f0, sN, fs, nK, t)
  135. # two freqs
  136. #h = minHarmonic2( f0+1e-2, f1-1e-2, sN, fs, nK, t)
  137. #h = harmonicEuler2( f0, f1, sN, fs, nK, t)
  138. plt.figure()
  139. plt.plot(t, sN, label="sN")
  140. #plt.plot(t, sN-h, label="sN-h")
  141. plt.plot(t, h, label='h')
  142. plt.title("harmonic")
  143. plt.legend()
  144. plt.figure()
  145. plt.plot(t, sN-sNc, label='true noise')
  146. plt.plot(t, h, label='harmonic removal')
  147. plt.plot(t, np.exp(-t/T2), label="nmr")
  148. plt.legend()
  149. plt.title("true noise")
  150. plt.show()