1
0
forked from IPKM/nmreval
2022-03-08 10:27:40 +01:00

112 lines
3.9 KiB
Python

from typing import Union
import numpy as np
from scipy.integrate import simps
from scipy.special import gammaln
from .base import Distribution
class GGAlpha(Distribution):
name = r'General \Gamma (\alpha-process)'
parameter = [r'\alpha', r'\beta']
@staticmethod
def correlation(t, tau0, *args):
logtau0 = np.log(tau0)
logtau = np.linspace(logtau0-20, logtau0+20, num=4001)
taus = np.exp(logtau)
gtau = GGAlpha.distribution(taus, tau0, *args)
ret_val = np.array([simps(np.exp(-xvals/taus)*gtau, logtau) for xvals in t])
return ret_val
@staticmethod
def susceptibility(omega, tau0, *args):
logtau0 = np.log(tau0)
logtau = np.linspace(logtau0 - 20, logtau0 + 20, num=4001)
taus = np.exp(logtau)
gtau = GGAlpha.distribution(taus, tau0, *args)
ret_val = np.array([simps(1/(1+1j*xvals*taus) * gtau, logtau) for xvals in omega])
return ret_val
@staticmethod
def distribution(taus: Union[float, np.ndarray], tau: float, *args) -> Union[float, np.ndarray]:
alpha, beta = args
b_to_a = beta / alpha
norm = np.exp(gammaln(b_to_a) - b_to_a * np.log(b_to_a)) / alpha
t_to_t0 = taus / tau
ret_val = np.exp(-b_to_a * t_to_t0 ** alpha) * t_to_t0 ** beta
return ret_val / norm
class GGAlphaEW(Distribution):
name = r'General \Gamma (\alpha-process + EW)'
parameter = [r'\alpha', r'\beta']
@staticmethod
def correlation(t, tau0, *args):
logtau0 = np.log(tau0)
logtau = np.linspace(logtau0-20, logtau0+20, num=4001)
taus = np.exp(logtau)
gtau = GGAlphaEW.distribution(taus, tau0, *args)
# wir integrieren ueber lntau, nicht tau
ret_val = np.array([simps(np.exp(-xvals/taus)*gtau, logtau) for xvals in t])
return ret_val
@staticmethod
def susceptibility(omega, tau0, *args):
logtau0 = np.log(tau0)
logtau = np.linspace(logtau0 - 20, logtau0 + 20, num=4001)
taus = np.exp(logtau)
gtau = GGAlphaEW.distribution(taus, tau0, *args)
ret_val = np.array([simps(1/(1+1j*xvals*taus) * gtau, logtau) for xvals in omega])
return ret_val
@staticmethod
def distribution(tau: Union[float, np.ndarray], tau0: float, *args) -> Union[float, np.ndarray]:
alpha, beta, sigma, gamma = args
if gamma == beta:
return GGAlpha.distribution(tau, tau0, alpha, beta)
b_to_a = beta / alpha
g_to_a = gamma / alpha
t_to_t0 = tau / tau0
norm = (np.exp(gammaln(b_to_a)) + sigma**(gamma-beta) *
np.exp(gammaln(g_to_a) + (b_to_a-g_to_a)*np.log(b_to_a))) / np.exp(b_to_a*np.log(b_to_a)) / alpha
ret_val = np.exp(-b_to_a * t_to_t0**alpha) * t_to_t0**beta * (1 + (t_to_t0*sigma)**(gamma-beta))
return ret_val / norm
class GGBeta(Distribution):
name = r'General \Gamma (\beta-process)'
parameter = [r'\alpha', r'\beta']
@staticmethod
def correlation(t, tau0, *args):
logtau0 = np.log(tau0)
logtau = np.linspace(logtau0-20, logtau0+20, num=4001)
taus = np.exp(logtau)
gtau = GGBeta.distribution(taus, tau0, *args)
# wir integrieren ueber lntau, nicht tau
ret_val = np.array([simps(np.exp(-xvals/taus)*gtau, logtau) for xvals in t])
return ret_val
@staticmethod
def susceptibility(omega, tau0, *args):
logtau0 = np.log(tau0)
logtau = np.linspace(logtau0 - 20, logtau0 + 20, num=4001)
taus = np.exp(logtau)
gtau = GGBeta.distribution(taus, tau0, *args)
ret_val = np.array([simps(1/(1+1j*xvals*taus) * gtau, logtau) for xvals in omega])
return ret_val
@staticmethod
def distribution(tau: Union[float, np.ndarray], tau0: float, *args) -> Union[float, np.ndarray]:
a, b = args
norm = a * (1+b) * np.sin(np.pi*b/(1+b)) * b**(b/(1+b)) / np.pi
ret_val = b * (tau/tau0)**a + (tau/tau0)**(-a*b)
return norm / ret_val