Package lib :: Package dispersion :: Module ns_cpmg_2site_expanded
[hide private]
[frames] | no frames]

Source Code for Module lib.dispersion.ns_cpmg_2site_expanded

  1  ############################################################################### 
  2  #                                                                             # 
  3  # Copyright (C) 2000-2001 Nikolai Skrynnikov                                  # 
  4  # Copyright (C) 2000-2001 Martin Tollinger                                    # 
  5  # Copyright (C) 2010-2013 Paul Schanda (https://web.archive.org/web/https://gna.org/users/pasa)           # 
  6  # Copyright (C) 2013 Mathilde Lescanne                                        # 
  7  # Copyright (C) 2013 Dominique Marion                                         # 
  8  # Copyright (C) 2013-2014 Edward d'Auvergne                                   # 
  9  #                                                                             # 
 10  # This file is part of the program relax (http://www.nmr-relax.com).          # 
 11  #                                                                             # 
 12  # This program is free software: you can redistribute it and/or modify        # 
 13  # it under the terms of the GNU General Public License as published by        # 
 14  # the Free Software Foundation, either version 3 of the License, or           # 
 15  # (at your option) any later version.                                         # 
 16  #                                                                             # 
 17  # This program is distributed in the hope that it will be useful,             # 
 18  # but WITHOUT ANY WARRANTY without even the implied warranty of               # 
 19  # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the               # 
 20  # GNU General Public License for more details.                                # 
 21  #                                                                             # 
 22  # You should have received a copy of the GNU General Public License           # 
 23  # along with this program.  If not, see <http://www.gnu.org/licenses/>.       # 
 24  #                                                                             # 
 25  ############################################################################### 
 26   
 27  # Module docstring. 
 28  """The numerical fit of 2-site Bloch-McConnell equations for CPMG-type experiments, the U{NS CPMG 2-site expanded<http://wiki.nmr-relax.com/NS_CPMG_2-site_expanded>} model. 
 29   
 30  Description 
 31  =========== 
 32   
 33  This function is exact, just as the explicit Bloch-McConnell numerical treatments.  It comes from a Maple derivation based on the Bloch-McConnell equations.  It is much faster than the numerical Bloch-McConnell solution.  It was derived by Nikolai Skrynnikov and is provided with his permission. 
 34   
 35   
 36  Code origin 
 37  =========== 
 38   
 39  The code originates as optimization function number 5 from the fitting_main_kex.py script from Mathilde Lescanne, Paul Schanda, and Dominique Marion (see U{http://thread.gmane.org/gmane.science.nmr.relax.devel/4138}, U{https://web.archive.org/web/https://gna.org/task/?7712#comment2} and U{https://web.archive.org/web/https://gna.org/support/download.php?file_id=18262}). 
 40   
 41  Links to the copyright licensing agreements from all authors are: 
 42   
 43      - Nikolai Skrynnikov, U{http://article.gmane.org/gmane.science.nmr.relax.devel/4279}, 
 44      - Martin Tollinger, U{http://article.gmane.org/gmane.science.nmr.relax.devel/4276}, 
 45      - Paul Schanda, U{http://article.gmane.org/gmane.science.nmr.relax.devel/4271}, 
 46      - Mathilde Lescanne, U{http://article.gmane.org/gmane.science.nmr.relax.devel/4138}, 
 47      - Dominique Marion, U{http://article.gmane.org/gmane.science.nmr.relax.devel/4157}. 
 48   
 49  The complex path of the code from the original Maple to relax can be described as: 
 50   
 51      - p3.analytical (Maple input text file at U{https://web.archive.org/web/https://gna.org/task/?7712#comment8}), 
 52      - Automatically generated FORTRAN, 
 53      - Manually converted to Matlab by Nikolai (sim_all.tar at U{https://web.archive.org/web/https://gna.org/task/?7712#comment5}) 
 54      - Manually converted to Python by Paul, Mathilde, and Dominique (fitting_main.py at U{https://web.archive.org/web/https://gna.org/task/?7712#comment1}) 
 55      - Converted into Python code for relax (here). 
 56   
 57  For reference, the original Maple script written by Nikolai for the expansion of the equations is:: 
 58   
 59      with(linalg): 
 60      with(tensor): 
 61      #Ka:=30; 
 62      #Kb:=1200; 
 63      #dW:=300; 
 64      #N:=2; 
 65      #tcp:=0.040/N; 
 66       
 67      Ksym:=sqrt(Ka*Kb); 
 68      #dX:=(Ka-Kb+I*dw)/2;        # Ra=Rb 
 69      dX:=((Ra-Rb)+(Ka-Kb)+I*dw)/2; 
 70       
 71      L:=([[-dX, Ksym], [Ksym, dX]]); 
 72        
 73      # in the end everything is multiplied by exp(-0.5*(Ra+Rb+Ka+Kb)*(Tc+2*tpalmer)) 
 74      # where 0.5*(Ra+Rb) is the same as Rinf, and (Ka+Kb) is kex. 
 75       
 76      y:=eigenvects(L); 
 77      TP1:=array([[y[1][3][1][1],y[2][3][1][1]],[y[1][3][1][2],y[2][3][1][2]]]); 
 78      iTP1:=inverse(TP1); 
 79      P1:=array([[exp(y[1][1]*tcp/2),0],[0,exp(y[2][1]*tcp/2)]]); 
 80       
 81      P1palmer:=array([[exp(y[1][1]*tpalmer),0],[0,exp(y[2][1]*tpalmer)]]); 
 82       
 83      TP2:=map(z->conj(z),TP1); 
 84      iTP2:=map(z->conj(z),iTP1); 
 85      P2:=array([[exp(conj(y[1][1])*tcp),0],[0,exp(conj(y[2][1])*tcp)]]); 
 86       
 87      P2palmer:=array([[exp(conj(y[1][1])*tpalmer),0],[0,exp(conj(y[2][1])*tpalmer)]]); 
 88       
 89      cP1:=evalm(TP1&*P1&*iTP1); 
 90      cP2:=evalm(TP2&*P2&*iTP2); 
 91       
 92      cP1palmer:=evalm(TP1&*P1palmer&*iTP1); 
 93      cP2palmer:=evalm(TP2&*P2palmer&*iTP2); 
 94       
 95      Ps:=evalm(cP1&*cP2&*cP1); 
 96      # Ps is symmetric; cf. simplify(Ps[1,2]-Ps[2,1]); 
 97      Pspalmer:=evalm(cP2palmer&*cP1palmer); 
 98       
 99       
100      dummy:=array([[a,b],[b,c]]); 
101      x:=eigenvects(dummy); 
102      TG1:=array([[x[1][3][1][1],x[2][3][1][1]],[x[1][3][1][2],x[2][3][1][2]]]); 
103      iTG1:=inverse(TG1); 
104      G1:=array([[x[1][1]^(N/4),0],[0,x[2][1]^(N/4)]]); 
105      GG1:=evalm(TG1&*G1&*iTG1); 
106      GG2:=map(z->conj(z),GG1); 
107       
108      cGG:=evalm(GG2&*Pspalmer&*GG1); 
109       
110      #s0:=array([Kb, Ka]); 
111      s0:=array([sqrt(Kb),sqrt(Ka)]);     # accounts for exchange symmetrization 
112      st:=evalm(cGG&*s0); 
113      #obs:=(1/(Ka+Kb))*st[1]; 
114      obs:=(sqrt(Kb)/(Ka+Kb))*st[1];  # accounts for exchange symmetrization 
115       
116      obs1:=eval(obs,[a=Ps[1,1],b=Ps[1,2],c=Ps[2,2]]); 
117      #obs2:=simplify(obs1): 
118       
119      print(obs1): 
120       
121      cGGref:=evalm(Pspalmer); 
122      stref:=evalm(cGGref&*s0); 
123      obsref:=(sqrt(Kb)/(Ka+Kb))*stref[1];  # accounts for exchange symmetrization 
124       
125      print(obsref): 
126       
127      writeto(result_test): 
128       
129      fortran([intensity=obs1, intensity_ref=obsref], optimized): 
130   
131   
132  Also for reference, the Matlab code from Nikolai and Martin manually converted from the automatically generated FORTRAN from the previous script into the funNikolai.m file is:: 
133   
134      function residual = funNikolai(optpar) 
135   
136      % extended Carver's equation derived via Maple, Ra-Rb = 0 by Skrynnikov 
137       
138      global nu_0 x y Rcalc rms nfields 
139      global Tc 
140       
141      Rcalc=zeros(nfields,size(x,2)); 
142       
143      tau_ex=optpar(1); 
144      pb=optpar(2); 
145       
146      pa=1-pb; 
147      kex=1/tau_ex; 
148      Ka=kex*pb; 
149      Kb=kex*pa; 
150       
151      nu_cpmg=x; 
152      tcp=1./(2*nu_cpmg); 
153      N=round(Tc./tcp); 
154       
155      for k=1:nfields 
156          dw=2*pi*nu_0(k)*optpar(3); 
157          Rinf=optpar(3+k); 
158           
159          t3 = i; 
160          t4 = t3*dw; 
161          t5 = Kb^2; 
162          t8 = 2*t3*Kb*dw; 
163          t10 = 2*Kb*Ka; 
164          t11 = dw^2; 
165          t14 = 2*t3*Ka*dw; 
166          t15 = Ka^2; 
167          t17 = sqrt(t5-t8+t10-t11+t14+t15); 
168          t21 = exp((-Kb+t4-Ka+t17)*tcp/4); 
169          t22 = 1/t17; 
170          t28 = exp((-Kb+t4-Ka-t17)*tcp/4); 
171          t31 = t21*t22*Ka-t28*t22*Ka; 
172          t33 = sqrt(t5+t8+t10-t11-t14+t15); 
173          t34 = Kb+t4-Ka+t33; 
174          t37 = exp((-Kb-t4-Ka+t33)*tcp/2); 
175          t39 = 1/t33; 
176          t41 = Kb+t4-Ka-t33; 
177          t44 = exp((-Kb-t4-Ka-t33)*tcp/2); 
178          t47 = t34*t37*t39/2-t41*t44*t39/2; 
179          t49 = Kb-t4-Ka-t17; 
180          t51 = t21*t49*t22; 
181          t52 = Kb-t4-Ka+t17; 
182          t54 = t28*t52*t22; 
183          t55 = -t51+t54; 
184          t60 = t37*t39*Ka-t44*t39*Ka; 
185          t62 = t31.*t47+t55.*t60/2; 
186          t63 = 1/Ka; 
187          t68 = -t52*t63*t51/2+t49*t63*t54/2; 
188          t69 = t62.*t68/2; 
189          t72 = t37*t41*t39; 
190          t76 = t44*t34*t39; 
191          t78 = -t34*t63*t72/2+t41*t63*t76/2; 
192          t80 = -t72+t76; 
193          t82 = t31.*t78/2+t55.*t80/4; 
194          t83 = t82.*t55/2; 
195          t88 = t52*t21*t22/2-t49*t28*t22/2; 
196          t91 = t88.*t47+t68.*t60/2; 
197          t92 = t91.*t88; 
198          t95 = t88.*t78/2+t68.*t80/4; 
199          t96 = t95.*t31; 
200          t97 = t69+t83; 
201          t98 = t97.^2; 
202          t99 = t92+t96; 
203          t102 = t99.^2; 
204          t108 = t62.*t88+t82.*t31; 
205          t112 = sqrt(t98-2*t99.*t97+t102+4*(t91.*t68/2+t95.*t55/2).*t108); 
206          t113 = t69+t83-t92-t96-t112; 
207          t115 = N/2; 
208          t116 = (t69/2+t83/2+t92/2+t96/2+t112/2).^t115; 
209          t118 = 1./t112; 
210          t120 = t69+t83-t92-t96+t112; 
211          t122 = (t69/2+t83/2+t92/2+t96/2-t112/2).^t115; 
212          t127 = 1./t108; 
213          t139 = 1/(Ka+Kb)*((-t113.*t116.*t118/2+t120.*t122.*t118/2)*Kb+(-t113.*t127.*t116.*t120.*t118/2+t120.*t127.*t122.*t113.*t118/2)*Ka/2); 
214           
215          intensity0 = pa;                             % pA 
216          intensity = real(t139)*exp(-Tc*Rinf);        % that's "homogeneous" relaxation 
217          Rcalc(k,:)=(1/Tc)*log(intensity0./intensity);  
218           
219      end 
220       
221      if (size(Rcalc)==size(y)) 
222          residual=sum(sum((y-Rcalc).^2));  
223          rms=sqrt(residual/(size(y,1)*size(y,2))); 
224      end 
225   
226   
227  Links 
228  ===== 
229   
230  More information on the NS CPMG 2-site expanded model can be found in the: 
231   
232      - U{relax wiki<http://wiki.nmr-relax.com/NS_CPMG_2-site_expanded>}, 
233      - U{relax manual<http://www.nmr-relax.com/manual/NS_2_site_expanded_CPMG_model.html>}, 
234      - U{relaxation dispersion page of the relax website<http://www.nmr-relax.com/analyses/relaxation_dispersion.html#NS_CPMG_2-site_expanded>}. 
235  """ 
236   
237  # Dependency check module. 
238  import dep_check 
239   
240  # Python module imports. 
241  from math import log 
242  from numpy import add, conj, dot, exp, power, real, sqrt 
243   
244  # relax module imports. 
245  from lib.float import isNaN 
246   
247   
248 -def r2eff_ns_cpmg_2site_expanded(r20=None, pA=None, dw=None, k_AB=None, k_BA=None, relax_time=None, inv_relax_time=None, tcp=None, back_calc=None, num_points=None, num_cpmg=None):
249 """The 2-site numerical solution to the Bloch-McConnell equation using complex conjugate matrices. 250 251 This function calculates and stores the R2eff values. 252 253 254 @keyword r20: The R2 value for both states A and B in the absence of exchange. 255 @type r20: float 256 @keyword pA: The population of state A. 257 @type pA: float 258 @keyword dw: The chemical exchange difference between states A and B in rad/s. 259 @type dw: float 260 @keyword k_AB: The rate of exchange from site A to B (rad/s). 261 @type k_AB: float 262 @keyword k_BA: The rate of exchange from site B to A (rad/s). 263 @type k_BA: float 264 @keyword relax_time: The total relaxation time period (in seconds). 265 @type relax_time: float 266 @keyword inv_relax_time: The inverse of the total relaxation time period (in inverse seconds). 267 @type inv_relax_time: float 268 @keyword tcp: The tau_CPMG times (1 / 4.nu1). 269 @type tcp: numpy rank-1 float array 270 @keyword back_calc: The array for holding the back calculated R2eff values. Each element corresponds to one of the CPMG nu1 frequencies. 271 @type back_calc: numpy rank-1 float array 272 @keyword num_points: The number of points on the dispersion curve, equal to the length of the tcp and back_calc arguments. 273 @type num_points: int 274 @keyword num_cpmg: The array of numbers of CPMG blocks. 275 @type num_cpmg: numpy int16, rank-1 array 276 """ 277 278 # Repeditive calculations. 279 half_tcp = 0.5 * tcp 280 k_AB_plus_k_BA = k_AB + k_BA 281 k_BA_minus_k_AB = k_BA - k_AB 282 283 # The expansion factors (in numpy array form for all dispersion points). 284 t3 = 1.j 285 t4 = t3 * dw 286 two_t4 = 2.0 * t4 287 t5 = k_BA**2 288 t8 = two_t4 * k_BA 289 t10 = 2.0 * k_BA * k_AB 290 t11 = dw**2 291 t14 = two_t4 * k_AB 292 t15 = k_AB**2 293 t5_t10_t11_t15 = t5 + t10 - t11 + t15 294 t8_t14 = t8 - t14 295 t17 = sqrt(t5_t10_t11_t15 - t8_t14) 296 297 k_AB_plus_k_BA_minus_t4 = k_AB_plus_k_BA - t4 298 t21 = exp((t17 - k_AB_plus_k_BA_minus_t4) * half_tcp) 299 t22 = 1.0/t17 300 t28 = exp(-(t17 + k_AB_plus_k_BA_minus_t4) * half_tcp) 301 t31 = t22*k_AB * (t21 - t28) 302 t33 = sqrt(t5_t10_t11_t15 + t8_t14) 303 304 k_AB_plus_k_BA_plus_t4 = k_AB_plus_k_BA + t4 305 k_BA_minus_k_AB_plus_t4 = k_BA_minus_k_AB + t4 306 t34 = k_BA_minus_k_AB_plus_t4 + t33 307 t37 = exp((t33 - k_AB_plus_k_BA_plus_t4) * tcp) 308 t39 = 1.0/t33 309 t41 = k_BA_minus_k_AB_plus_t4 - t33 310 t44 = exp(-(t33 + k_AB_plus_k_BA_plus_t4) * tcp) 311 t47 = 0.5*t39 * (t34*t37 - t41*t44) 312 313 k_BA_minus_k_AB_minus_t4 = k_BA_minus_k_AB - t4 314 t49 = k_BA_minus_k_AB_minus_t4 - t17 315 t51 = t21 * t49 * t22 316 t52 = k_BA_minus_k_AB_minus_t4 + t17 317 t54 = t28 * t52 * t22 318 t55 = -t51 + t54 319 t60 = 0.5*t39*k_AB * (t37 - t44) 320 t62 = t31*t47 + t55*t60 321 t63 = 1.0/k_AB 322 t68 = 0.5*t63 * (t49*t54 - t52*t51) 323 t69 = 0.5*t62 * t68 324 t72 = t37 * t41 * t39 325 t76 = t44 * t34 * t39 326 t78 = 0.5*t63 * (t41*t76 - t34*t72) 327 t80 = 0.5 * (t76 - t72) 328 t82 = 0.5 * (t31*t78 + t55*t80) 329 t83 = t82 * t55/2.0 330 t88 = 0.5 * t22 * (t52*t21 - t49*t28) 331 t91 = t88 * t47 + t68*t60 332 t92 = t91 * t88 333 t95 = 0.5 * (t88*t78 + t68*t80) 334 t96 = t95 * t31 335 t97 = t69 + t83 336 t98 = t97**2 337 t99 = t92 + t96 338 t102 = t99**2 339 t108 = t62 * t88 + t82 * t31 340 t112 = sqrt(t98 - 2.0 * t99 * t97 + t102 + 2.0 * (t91 * t68 + t95 * t55) * t108) 341 t97_t99 = t97 + t99 342 t97_nt99 = t97 - t99 343 t113 = t97_nt99 - t112 344 t115 = num_cpmg 345 t116 = power(0.5*(t97_t99 + t112), t115) 346 t118 = 1.0/t112 347 t120 = t97_nt99 + t112 348 t122 = power(0.5*(t97_t99 - t112), t115) 349 t127 = 0.5/t108 350 t120_t122 = t120*t122 351 t139 = 0.5/(k_AB + k_BA) * ((t120_t122 - t113*t116)*t118*k_BA + (t120_t122 - t116*t120)*t127*t113*t118*k_AB) 352 353 # Calculate the initial and final peak intensities. 354 intensity0 = pA 355 intensity = t139.real * exp(-relax_time * r20) 356 357 # The magnetisation vector. 358 Mx = intensity / intensity0 359 360 # Calculate the R2eff using a two-point approximation, i.e. assuming that the decay is mono-exponential, and store it for each dispersion point. 361 for i in range(num_points): 362 if Mx[i] <= 0.0 or isNaN(Mx[i]): 363 back_calc[i] = 1e99 364 else: 365 back_calc[i]= -inv_relax_time * log(Mx[i])
366