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.

1206 lines
47 KiB

1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
1 year ago
  1. import numpy as np
  2. from uncertainties import ufloat
  3. import lmfit
  4. from lmfit.models import (ConstantModel, ComplexConstantModel, LinearModel, QuadraticModel,
  5. PolynomialModel, SineModel, GaussianModel, Gaussian2dModel, LorentzianModel,
  6. SplitLorentzianModel, VoigtModel, PseudoVoigtModel, MoffatModel,
  7. Pearson7Model, StudentsTModel, BreitWignerModel, LognormalModel,
  8. DampedOscillatorModel, ExponentialGaussianModel, SkewedGaussianModel,
  9. SkewedVoigtModel, ThermalDistributionModel, DoniachModel, PowerLawModel,
  10. ExponentialModel, StepModel, RectangleModel, ExpressionModel, DampedHarmonicOscillatorModel)
  11. from lmfit.models import (guess_from_peak, guess_from_peak2d, fwhm_expr, height_expr,
  12. update_param_vals)
  13. from lmfit.lineshapes import (not_zero, breit_wigner, damped_oscillator, dho, doniach,
  14. expgaussian, exponential, gaussian, gaussian2d,
  15. linear, lognormal, lorentzian, moffat, parabolic,
  16. pearson7, powerlaw, pvoigt, rectangle, sine,
  17. skewed_gaussian, skewed_voigt, split_lorentzian, step,
  18. students_t, thermal_distribution, tiny, voigt)
  19. from lmfit import Model
  20. import numpy as np
  21. from numpy import (arctan, copysign, cos, exp, isclose, isnan, log, pi, real,
  22. sin, sqrt, where)
  23. from scipy.special import erf, erfc
  24. from scipy.special import gamma as gamfcn
  25. from scipy.special import wofz
  26. from scipy.optimize import curve_fit
  27. import xarray as xr
  28. import copy
  29. log2 = log(2)
  30. s2pi = sqrt(2*pi)
  31. s2 = sqrt(2.0)
  32. def gaussianWithOffset(x, amplitude=1.0, center=0.0, sigma=1.0, offset=0.0):
  33. """Return a 1-dimensional Gaussian function with an offset.
  34. gaussian(x, amplitude, center, sigma) =
  35. (amplitude/(s2pi*sigma)) * exp(-(1.0*x-center)**2 / (2*sigma**2))
  36. """
  37. return ((amplitude/(max(tiny, s2pi*sigma)))
  38. * exp(-(1.0*x-center)**2 / max(tiny, (2*sigma**2))) + offset)
  39. def lorentzianWithOffset(x, amplitude=1.0, center=0.0, sigma=1.0, offset=0.0):
  40. return ((amplitude/(1 + ((1.0*x-center)/max(tiny, sigma))**2))
  41. / max(tiny, (pi*sigma)) + offset)
  42. def exponentialWithOffset(x, amplitude=1.0, decay=1.0, offset=0.0):
  43. decay = not_zero(decay)
  44. return amplitude * exp(-x/decay) + offset
  45. def expansion(x, amplitude=1.0, offset=0.0):
  46. return np.sqrt(amplitude*x*x + offset)
  47. def dampingOscillation(x, center=0, amplitude=1.0, frequency=1.0, decay=1.0, offset=0.0):
  48. return amplitude * np.exp(-decay*x)*np.sin(2*np.pi*frequency*(x-center)) + offset
  49. def two_gaussian2d(x, y=0.0, A_amplitude=1.0, A_centerx=0.0, A_centery=0.0, A_sigmax=1.0,
  50. A_sigmay=1.0, B_amplitude=1.0, B_centerx=0.0, B_centery=0.0, B_sigmax=1.0,
  51. B_sigmay=1.0):
  52. """Return a 2-dimensional Gaussian function.
  53. gaussian2d(x, y, amplitude, centerx, centery, sigmax, sigmay) =
  54. amplitude/(2*pi*sigmax*sigmay) * exp(-(x-centerx)**2/(2*sigmax**2)
  55. -(y-centery)**2/(2*sigmay**2))
  56. """
  57. z = A_amplitude*(gaussian(x, amplitude=1, center=A_centerx, sigma=A_sigmax) *
  58. gaussian(y, amplitude=1, center=A_centery, sigma=A_sigmay))
  59. z += B_amplitude*(gaussian(x, amplitude=1, center=B_centerx, sigma=B_sigmax) *
  60. gaussian(y, amplitude=1, center=B_centery, sigma=B_sigmay))
  61. return z
  62. def ThomasFermi_2d(x, y=0.0, centerx=0.0, centery=0.0, amplitude=1.0, sigmax=1.0, sigmay=1.0):
  63. res = (1- ((x-centerx)/(sigmax))**2 - ((y-centery)/(sigmay))**2)**(3 / 2)
  64. return amplitude * 5 / 2 / np.pi / max(tiny, sigmax * sigmay) * np.where(res > 0, res, 0)
  65. def polylog(power, numerator):
  66. order = 20
  67. dataShape = numerator.shape
  68. numerator = np.tile(numerator, (order, 1))
  69. numerator = np.power(numerator.T, np.arange(1, order+1)).T
  70. denominator = np.arange(1, order+1)
  71. denominator = np.tile(denominator, (dataShape[0], 1))
  72. denominator = denominator.T
  73. data = numerator/ np.power(denominator, power)
  74. return np.sum(data, axis=0)
  75. def polylog2_2d(x, y=0.0, centerx=0.0, centery=0.0, amplitude=1.0, sigmax=1.0, sigmay=1.0):
  76. ## Approximation of the polylog function with 2D gaussian as argument. -> discribes the thermal part of the cloud
  77. return amplitude / 2 / np.pi / 1.20206 / max(tiny, sigmax * sigmay) * polylog(2, np.exp( -((x-centerx)**2/(2 * (sigmax)**2))-((y-centery)**2/( 2 * (sigmay)**2)) ))
  78. def density_profile_BEC_2d(x, y=0.0, BEC_amplitude=1.0, thermal_amplitude=1.0, BEC_centerx=0.0, BEC_centery=0.0, thermal_centerx=0.0, thermal_centery=0.0,
  79. BEC_sigmax=1.0, BEC_sigmay=1.0, thermal_sigmax=1.0, thermal_sigmay=1.0):
  80. return ThomasFermi_2d(x=x, y=y, centerx=BEC_centerx, centery=BEC_centery,
  81. amplitude=BEC_amplitude, sigmax=BEC_sigmax, sigmay=BEC_sigmay
  82. ) + polylog2_2d(x=x, y=y, centerx=thermal_centerx, centery=thermal_centery,
  83. amplitude=thermal_amplitude, sigmax=thermal_sigmax, sigmay=thermal_sigmay)
  84. class GaussianWithOffsetModel(Model):
  85. fwhm_factor = 2*np.sqrt(2*np.log(2))
  86. height_factor = 1./np.sqrt(2*np.pi)
  87. def __init__(self, independent_vars=['x'], nan_policy='raise', prefix='', name=None, **kwargs):
  88. kwargs.update({'prefix': prefix, 'nan_policy': nan_policy,
  89. 'independent_vars': independent_vars})
  90. super().__init__(gaussianWithOffset, **kwargs)
  91. self._set_paramhints_prefix()
  92. def _set_paramhints_prefix(self):
  93. self.set_param_hint('sigma', min=0)
  94. self.set_param_hint('fwhm', expr=fwhm_expr(self))
  95. self.set_param_hint('height', expr=height_expr(self))
  96. def guess(self, data, x, negative=False, **kwargs):
  97. offset = np.min(data)
  98. data = data - offset
  99. pars = guess_from_peak(self, data, x, negative)
  100. pars.add('offset', value=offset)
  101. return update_param_vals(pars, self.prefix, **kwargs)
  102. class LorentzianWithOffsetModel(Model):
  103. fwhm_factor = 2.0
  104. height_factor = 1./np.pi
  105. def __init__(self, independent_vars=['x'], prefix='', nan_policy='raise',
  106. **kwargs):
  107. kwargs.update({'prefix': prefix, 'nan_policy': nan_policy,
  108. 'independent_vars': independent_vars})
  109. super().__init__(lorentzianWithOffset, **kwargs)
  110. self._set_paramhints_prefix()
  111. def _set_paramhints_prefix(self):
  112. self.set_param_hint('sigma', min=0)
  113. self.set_param_hint('fwhm', expr=fwhm_expr(self))
  114. self.set_param_hint('height', expr=height_expr(self))
  115. def guess(self, data, x, negative=False, **kwargs):
  116. """Estimate initial model parameter values from data."""
  117. offset = np.min(data)
  118. data = data - offset
  119. pars = guess_from_peak(self, data, x, negative, ampscale=1.25)
  120. pars.add('offset', value=offset)
  121. return update_param_vals(pars, self.prefix, **kwargs)
  122. class ExponentialWithOffsetModel(Model):
  123. def __init__(self, independent_vars=['x'], prefix='', nan_policy='raise',
  124. **kwargs):
  125. kwargs.update({'prefix': prefix, 'nan_policy': nan_policy,
  126. 'independent_vars': independent_vars})
  127. super().__init__(exponentialWithOffset, **kwargs)
  128. def guess(self, data, x, **kwargs):
  129. """Estimate initial model parameter values from data."""
  130. offset = np.min(data)
  131. data = data - offset
  132. try:
  133. sval, oval = np.polyfit(x, np.log(abs(data)+1.e-15), 1)
  134. except TypeError:
  135. sval, oval = 1., np.log(abs(max(data)+1.e-9))
  136. pars = self.make_params(amplitude=np.exp(oval), decay=-1.0/sval)
  137. pars.add('offset', value=offset)
  138. return update_param_vals(pars, self.prefix, **kwargs)
  139. class ExpansionModel(Model):
  140. def __init__(self, independent_vars=['x'], prefix='', nan_policy='raise',
  141. **kwargs):
  142. kwargs.update({'prefix': prefix, 'nan_policy': nan_policy,
  143. 'independent_vars': independent_vars})
  144. super().__init__(expansion, **kwargs)
  145. def guess(self, data, x, **kwargs):
  146. """Estimate initial model parameter values from data."""
  147. popt1, pcov1 = curve_fit(expansion, x, data)
  148. pars = self.make_params(amplitude=popt1[0], offset=popt1[1])
  149. return update_param_vals(pars, self.prefix, **kwargs)
  150. class DampingOscillationModel(Model):
  151. def __init__(self, independent_vars=['x'], prefix='', nan_policy='raise',
  152. **kwargs):
  153. kwargs.update({'prefix': prefix, 'nan_policy': nan_policy,
  154. 'independent_vars': independent_vars})
  155. super().__init__(dampingOscillation, **kwargs)
  156. def guess(self, data, x, **kwargs):
  157. """Estimate initial model parameter values from data."""
  158. try:
  159. popt1, pcov1 = curve_fit(dampingOscillation, x, data, np.array(0, 5, 5e2, 1e3, 16))
  160. pars = self.make_params(center=popt1[0], amplitude=popt1[1], frequency=popt1[2], decay=popt1[3], offset=popt1[4])
  161. except:
  162. pars = self.make_params(center=0, amplitude=5.0, frequency=5e2, decay=1.0e3, offset=16.0)
  163. return update_param_vals(pars, self.prefix, **kwargs)
  164. class TwoGaussian2dModel(Model):
  165. fwhm_factor = 2*np.sqrt(2*np.log(2))
  166. height_factor = 1./2*np.pi
  167. def __init__(self, independent_vars=['x', 'y'], prefix='', nan_policy='raise',
  168. **kwargs):
  169. kwargs.update({'prefix': prefix, 'nan_policy': nan_policy,
  170. 'independent_vars': independent_vars})
  171. self.helperModel = Gaussian2dModel()
  172. super().__init__(two_gaussian2d, **kwargs)
  173. self._set_paramhints_prefix()
  174. def _set_paramhints_prefix(self):
  175. self.set_param_hint('delta', value=-1, max=0)
  176. self.set_param_hint('A_sigmax', expr=f'{self.prefix}delta + {self.prefix}B_sigmax')
  177. def guess(self, data, x, y, negative=False, **kwargs):
  178. pars_guess = guess_from_peak2d(self.helperModel, data, x, y, negative)
  179. pars = self.make_params(A_amplitude=pars_guess['amplitude'].value, A_centerx=pars_guess['centerx'].value, A_centery=pars_guess['centery'].value,
  180. A_sigmax=pars_guess['sigmax'].value, A_sigmay=pars_guess['sigmay'].value,
  181. B_amplitude=pars_guess['amplitude'].value, B_centerx=pars_guess['centerx'].value, B_centery=pars_guess['centery'].value,
  182. B_sigmax=pars_guess['sigmax'].value, B_sigmay=pars_guess['sigmay'].value)
  183. pars[f'{self.prefix}A_sigmax'].set(expr=f'delta + {self.prefix}B_sigmax')
  184. pars.add(f'{self.prefix}delta', value=-1, max=0, min=-np.inf, vary=True)
  185. pars[f'{self.prefix}A_sigmay'].set(min=0.0)
  186. pars[f'{self.prefix}B_sigmax'].set(min=0.0)
  187. pars[f'{self.prefix}B_sigmay'].set(min=0.0)
  188. return pars
  189. class Polylog22dModel(Model):
  190. fwhm_factor = 2*np.sqrt(2*np.log(2))
  191. height_factor = 1./2*np.pi
  192. def __init__(self, independent_vars=['x', 'y'], prefix='', nan_policy='raise',
  193. **kwargs):
  194. kwargs.update({'prefix': prefix, 'nan_policy': nan_policy,
  195. 'independent_vars': independent_vars})
  196. super().__init__(polylog2_2d, **kwargs)
  197. self._set_paramhints_prefix()
  198. def _set_paramhints_prefix(self):
  199. self.set_param_hint('Rx', min=0)
  200. self.set_param_hint('Ry', min=0)
  201. def guess(self, data, x, y, negative=False, **kwargs):
  202. """Estimate initial model parameter values from data."""
  203. pars = guess_from_peak2d(self, data, x, y, negative)
  204. return update_param_vals(pars, self.prefix, **kwargs)
  205. class ThomasFermi2dModel(Model):
  206. fwhm_factor = 1
  207. height_factor = 0.5
  208. def __init__(self, independent_vars=['x', 'y'], prefix='', nan_policy='raise',
  209. **kwargs):
  210. kwargs.update({'prefix': prefix, 'nan_policy': nan_policy,
  211. 'independent_vars': independent_vars})
  212. super().__init__(ThomasFermi_2d, **kwargs)
  213. self._set_paramhints_prefix()
  214. def _set_paramhints_prefix(self):
  215. self.set_param_hint('Rx', min=0)
  216. self.set_param_hint('Ry', min=0)
  217. def guess(self, data, x, y, negative=False, **kwargs):
  218. """Estimate initial model parameter values from data."""
  219. pars = guess_from_peak2d(self, data, x, y, negative)
  220. # amplitude = pars['amplitude'].value
  221. # simgax = pars['sigmax'].value
  222. # sigmay = pars['sigmay'].value
  223. # pars['amplitude'].set(value=amplitude/s2pi/simgax/sigmay)
  224. simgax = pars['sigmax'].value
  225. sigmay = pars['sigmay'].value
  226. pars['simgax'].set(value=simgax / 2.355)
  227. pars['sigmay'].set(value=sigmay / 2.355)
  228. return update_param_vals(pars, self.prefix, **kwargs)
  229. class DensityProfileBEC2dModel(Model):
  230. fwhm_factor = 2*np.sqrt(2*np.log(2))
  231. height_factor = 1./2*np.pi
  232. def __init__(self, independent_vars=['x', 'y'], prefix='', nan_policy='raise',
  233. **kwargs):
  234. kwargs.update({'prefix': prefix, 'nan_policy': nan_policy,
  235. 'independent_vars': independent_vars})
  236. super().__init__(density_profile_BEC_2d, **kwargs)
  237. self._set_paramhints_prefix()
  238. def _set_paramhints_prefix(self):
  239. # self.set_param_hint('BEC_sigmax', min=0)
  240. self.set_param_hint('deltax', min=0)
  241. self.set_param_hint('BEC_sigmax', expr=f'3 * {self.prefix}thermal_sigmax - {self.prefix}deltax')
  242. self.set_param_hint('BEC_sigmay', min=0)
  243. self.set_param_hint('thermal_sigmax', min=0)
  244. # self.set_param_hint('thermal_sigmay', min=0)
  245. self.set_param_hint('BEC_amplitude', min=0)
  246. self.set_param_hint('thermal_amplitude', min=0)
  247. self.set_param_hint('thermalAspectRatio', min=0.8, max=1.2)
  248. self.set_param_hint('thermal_sigmay', expr=f'{self.prefix}thermalAspectRatio * {self.prefix}thermal_sigmax')
  249. # self.set_param_hint('betax', value=0)
  250. # self.set_param_hint('BEC_centerx', expr=f'{self.prefix}thermal_sigmax - {self.prefix}betax')
  251. self.set_param_hint('condensate_fraction', expr=f'{self.prefix}BEC_amplitude / ({self.prefix}BEC_amplitude + {self.prefix}thermal_amplitude)')
  252. def guess(self, data, x, y, negative=False, pureBECThreshold=0.5, noBECThThreshold=0.0, **kwargs):
  253. """Estimate initial model parameter values from data."""
  254. fitModel = TwoGaussian2dModel()
  255. pars = fitModel.guess(data, x=x, y=y, negative=negative)
  256. pars['A_amplitude'].set(min=0)
  257. pars['B_amplitude'].set(min=0)
  258. pars['A_centerx'].set(min=pars['A_centerx'].value - 3 * pars['A_sigmax'],
  259. max=pars['A_centerx'].value + 3 * pars['A_sigmax'],)
  260. pars['A_centery'].set(min=pars['A_centery'].value - 3 * pars['A_sigmay'],
  261. max=pars['A_centery'].value + 3 * pars['A_sigmay'],)
  262. pars['B_centerx'].set(min=pars['B_centerx'].value - 3 * pars['B_sigmax'],
  263. max=pars['B_centerx'].value + 3 * pars['B_sigmax'],)
  264. pars['B_centery'].set(min=pars['B_centery'].value - 3 * pars['B_sigmay'],
  265. max=pars['B_centery'].value + 3 * pars['B_sigmay'],)
  266. fitResult = fitModel.fit(data, x=x, y=y, params=pars, **kwargs)
  267. pars_guess = fitResult.params
  268. BEC_amplitude = pars_guess['A_amplitude'].value
  269. thermal_amplitude = pars_guess['B_amplitude'].value
  270. pars = self.make_params(BEC_amplitude=BEC_amplitude,
  271. thermal_amplitude=thermal_amplitude,
  272. BEC_centerx=pars_guess['A_centerx'].value, BEC_centery=pars_guess['A_centery'].value,
  273. # BEC_sigmax=(pars_guess['A_sigmax'].value / 2.355),
  274. deltax = 3 * (pars_guess['B_sigmax'].value * s2) - (pars_guess['A_sigmax'].value / 2.355),
  275. BEC_sigmay=(pars_guess['A_sigmay'].value / 2.355),
  276. thermal_centerx=pars_guess['B_centerx'].value, thermal_centery=pars_guess['B_centery'].value,
  277. thermal_sigmax=(pars_guess['B_sigmax'].value * s2),
  278. thermalAspectRatio=(pars_guess['B_sigmax'].value * s2) / (pars_guess['B_sigmay'].value * s2)
  279. # thermal_sigmay=(pars_guess['B_sigmay'].value * s2)
  280. )
  281. nBEC = pars[f'{self.prefix}BEC_amplitude'] / 2 / np.pi / 5.546 / pars[f'{self.prefix}BEC_sigmay'] / pars[f'{self.prefix}BEC_sigmax']
  282. if (pars[f'{self.prefix}condensate_fraction']>0.95) and (np.max(data) > 1.05 * nBEC):
  283. temp = ((np.max(data) - nBEC) * s2pi * pars[f'{self.prefix}thermal_sigmay'] / pars[f'{self.prefix}thermal_sigmax'])
  284. if temp > pars[f'{self.prefix}BEC_amplitude']:
  285. pars[f'{self.prefix}thermal_amplitude'].set(value=pars[f'{self.prefix}BEC_amplitude'] / 2)
  286. else:
  287. pars[f'{self.prefix}thermal_amplitude'].set(value=temp * 10)
  288. if BEC_amplitude / (thermal_amplitude + BEC_amplitude) > pureBECThreshold:
  289. pars[f'{self.prefix}thermal_amplitude'].set(value=0)
  290. pars[f'{self.prefix}BEC_amplitude'].set(value=(thermal_amplitude + BEC_amplitude))
  291. if BEC_amplitude / (thermal_amplitude + BEC_amplitude) < noBECThThreshold:
  292. pars[f'{self.prefix}BEC_amplitude'].set(value=0)
  293. pars[f'{self.prefix}thermal_amplitude'].set(value=(thermal_amplitude + BEC_amplitude))
  294. pars[f'{self.prefix}BEC_centerx'].set(
  295. min=pars[f'{self.prefix}BEC_centerx'].value - 10 * pars[f'{self.prefix}BEC_sigmax'].value,
  296. max=pars[f'{self.prefix}BEC_centerx'].value + 10 * pars[f'{self.prefix}BEC_sigmax'].value,
  297. )
  298. pars[f'{self.prefix}thermal_centerx'].set(
  299. min=pars[f'{self.prefix}thermal_centerx'].value - 3 * pars[f'{self.prefix}thermal_sigmax'].value,
  300. max=pars[f'{self.prefix}thermal_centerx'].value + 3 * pars[f'{self.prefix}thermal_sigmax'].value,
  301. )
  302. pars[f'{self.prefix}BEC_centery'].set(
  303. min=pars[f'{self.prefix}BEC_centery'].value - 10 * pars[f'{self.prefix}BEC_sigmay'].value,
  304. max=pars[f'{self.prefix}BEC_centery'].value + 10 * pars[f'{self.prefix}BEC_sigmay'].value,
  305. )
  306. pars[f'{self.prefix}thermal_centery'].set(
  307. min=pars[f'{self.prefix}thermal_centery'].value - 3 * pars[f'{self.prefix}thermal_sigmay'].value,
  308. max=pars[f'{self.prefix}thermal_centery'].value + 3 * pars[f'{self.prefix}thermal_sigmay'].value,
  309. )
  310. pars[f'{self.prefix}BEC_sigmay'].set(
  311. max=5 * pars[f'{self.prefix}BEC_sigmay'].value,
  312. )
  313. pars[f'{self.prefix}thermal_sigmax'].set(
  314. max=5 * pars[f'{self.prefix}thermal_sigmax'].value,
  315. )
  316. return update_param_vals(pars, self.prefix, **kwargs)
  317. class NewFitModel(Model):
  318. def __init__(self, func, independent_vars=['x'], prefix='', nan_policy='raise',
  319. **kwargs):
  320. kwargs.update({'prefix': prefix, 'nan_policy': nan_policy,
  321. 'independent_vars': independent_vars})
  322. super().__init__(func, **kwargs)
  323. def guess(self, *args, **kwargs):
  324. return self.make_params()
  325. lmfit_models = {'Constant': ConstantModel,
  326. 'Complex Constant': ComplexConstantModel,
  327. 'Linear': LinearModel,
  328. 'Quadratic': QuadraticModel,
  329. 'Polynomial': PolynomialModel,
  330. 'Gaussian': GaussianModel,
  331. 'Gaussian-2D': Gaussian2dModel,
  332. 'Lorentzian': LorentzianModel,
  333. 'Split-Lorentzian': SplitLorentzianModel,
  334. 'Voigt': VoigtModel,
  335. 'PseudoVoigt': PseudoVoigtModel,
  336. 'Moffat': MoffatModel,
  337. 'Pearson7': Pearson7Model,
  338. 'StudentsT': StudentsTModel,
  339. 'Breit-Wigner': BreitWignerModel,
  340. 'Log-Normal': LognormalModel,
  341. 'Damped Oscillator': DampedOscillatorModel,
  342. 'Damped Harmonic Oscillator': DampedHarmonicOscillatorModel,
  343. 'Exponential Gaussian': ExponentialGaussianModel,
  344. 'Skewed Gaussian': SkewedGaussianModel,
  345. 'Skewed Voigt': SkewedVoigtModel,
  346. 'Thermal Distribution': ThermalDistributionModel,
  347. 'Doniach': DoniachModel,
  348. 'Power Law': PowerLawModel,
  349. 'Exponential': ExponentialModel,
  350. 'Step': StepModel,
  351. 'Rectangle': RectangleModel,
  352. 'Expression': ExpressionModel,
  353. 'Gaussian With Offset':GaussianWithOffsetModel,
  354. 'Lorentzian With Offset':LorentzianWithOffsetModel,
  355. 'Expansion':ExpansionModel,
  356. 'Damping Oscillation Model':DampingOscillationModel,
  357. 'Two Gaussian-2D':TwoGaussian2dModel,
  358. }
  359. class FitAnalyser():
  360. """This is a class integrated all the functions to do a fit.
  361. """
  362. def __init__(self, fitModel, fitDim=1, **kwargs) -> None:
  363. """Initialize function
  364. :param fitModel: The fitting model of fit function
  365. :type fitModel: lmfit Model
  366. :param fitDim: The dimension of the fit, defaults to 1
  367. :type fitDim: int, optional
  368. """
  369. if isinstance(fitModel, str):
  370. self.fitModel = lmfit_models[fitModel](**kwargs)
  371. else:
  372. self.fitModel = fitModel
  373. self.fitDim = fitDim
  374. def print_params_set_template(self, params=None):
  375. """Print a template to manually set the initial parameters of the fit
  376. :param params: An object of Parameters class to print, defaults to None
  377. :type params: lmfit Paraemters, optional
  378. """
  379. if params is None:
  380. params = self.fitModel.make_params()
  381. for key in params:
  382. res = "params.add("
  383. res += "name=\"" + key + "\", "
  384. if not params[key].expr is None:
  385. res += "expr=\"" + params[key].expr +"\")"
  386. else:
  387. res += "value=" + f'{params[key].value:3g}' + ", "
  388. if str(params[key].max)=="inf":
  389. res += "max=np.inf, "
  390. else:
  391. res += "max=" + f'{params[key].max:3g}' + ", "
  392. if str(params[key].min)=="-inf":
  393. res += "min=-np.inf, "
  394. else:
  395. res += "min=" + f'{params[key].min:3g}' + ", "
  396. res += "vary=" + str(params[key].vary) + ")"
  397. print(res)
  398. def _guess_1D(self, data, x, **kwargs):
  399. """Call the guess function of the 1D fit model to guess the initial value.
  400. :param data: The data to fit
  401. :type data: 1D numpy array
  402. :param x: The data of x axis
  403. :type x: 1D numpy array
  404. :return: The guessed initial parameters for the fit
  405. :rtype: lmfit Parameters
  406. """
  407. return self.fitModel.guess(data=data, x=x, **kwargs)
  408. def _guess_2D(self, data, x, y, **kwargs):
  409. """Call the guess function of the 2D fit model to guess the initial value.
  410. :param data: The flattened data to fit
  411. :type data: 1D numpy array
  412. :param x: The flattened data of x axis
  413. :type x: 1D numpy array
  414. :param y: The flattened data of y axis
  415. :type y: 1D numpy array
  416. :return: The guessed initial parameters for the fit
  417. :rtype: lmfit Parameters
  418. """
  419. data = data.flatten(order='F')
  420. return self.fitModel.guess(data=data, x=x, y=y, **kwargs)
  421. def guess(self, dataArray, x=None, y=None, guess_kwargs={}, input_core_dims=None, dask='parallelized', vectorize=True, keep_attrs=True, daskKwargs=None, **kwargs):
  422. """Call the guess function of the fit model to guess the initial value.
  423. :param dataArray: The data for the fit
  424. :type dataArray: xarray DataArray
  425. :param x: The name of x axis in data or the data of x axis, defaults to None
  426. :type x: str or numpy array, optional
  427. :param y: The name of y axis in data or the data of y axis, defaults to None
  428. :type y: str or numpy array, optional
  429. :param guess_kwargs: the keyworded arguments to send to the guess function, defaults to {}
  430. :type guess_kwargs: dict, optional
  431. :param input_core_dims: over write of the same argument in xarray.apply_ufunc, defaults to None
  432. :type input_core_dims: list or array like, optional
  433. :param dask: over write of the same argument in xarray.apply_ufunc,, defaults to 'parallelized'
  434. :type dask: str, optional
  435. :param vectorize: over write of the same argument in xarray.apply_ufunc, defaults to True
  436. :type vectorize: bool, optional
  437. :param keep_attrs: over write of the same argument in xarray.apply_ufunc, defaults to True
  438. :type keep_attrs: bool, optional
  439. :param daskKwargs: over write of the same argument in xarray.apply_ufunc, defaults to None
  440. :type daskKwargs: dict, optional
  441. :return: The guessed initial parameters for the fit
  442. :rtype: xarray DataArray
  443. """
  444. kwargs.update(
  445. {
  446. "dask": dask,
  447. "vectorize": vectorize,
  448. "input_core_dims": input_core_dims,
  449. 'keep_attrs': keep_attrs,
  450. }
  451. )
  452. if not daskKwargs is None:
  453. kwargs.update({"dask_gufunc_kwargs": daskKwargs})
  454. if input_core_dims is None:
  455. kwargs.update(
  456. {
  457. "input_core_dims": [['x']],
  458. }
  459. )
  460. if x is None:
  461. if 'x' in dataArray.dims:
  462. x = dataArray['x'].to_numpy()
  463. else:
  464. if isinstance(x, str):
  465. if input_core_dims is None:
  466. kwargs.update(
  467. {
  468. "input_core_dims": [[x]],
  469. }
  470. )
  471. x = dataArray[x].to_numpy()
  472. if self.fitDim == 1:
  473. guess_kwargs.update(
  474. {
  475. 'x':x,
  476. }
  477. )
  478. return xr.apply_ufunc(self._guess_1D, dataArray, kwargs=guess_kwargs,
  479. output_dtypes=[type(self.fitModel.make_params())],
  480. **kwargs
  481. )
  482. if self.fitDim == 2:
  483. if y is None:
  484. if 'y' in dataArray.dims:
  485. y = dataArray['y'].to_numpy()
  486. if input_core_dims is None:
  487. kwargs.update(
  488. {
  489. "input_core_dims": [['x', 'y']],
  490. }
  491. )
  492. else:
  493. if isinstance(y, str):
  494. kwargs["input_core_dims"][0] = np.append(kwargs["input_core_dims"][0], y)
  495. y = dataArray[y].to_numpy()
  496. elif input_core_dims is None:
  497. kwargs.update(
  498. {
  499. "input_core_dims": [['x', 'y']],
  500. }
  501. )
  502. _x, _y = np.meshgrid(x, y)
  503. _x = _x.flatten()
  504. _y = _y.flatten()
  505. # dataArray = dataArray.stack(_z=(kwargs["input_core_dims"][0][0], kwargs["input_core_dims"][0][1]))
  506. # kwargs["input_core_dims"][0] = ['_z']
  507. guess_kwargs.update(
  508. {
  509. 'x':_x,
  510. 'y':_y,
  511. }
  512. )
  513. return xr.apply_ufunc(self._guess_2D, dataArray, kwargs=guess_kwargs,
  514. output_dtypes=[type(self.fitModel.make_params())],
  515. **kwargs
  516. )
  517. def _fit_1D(self, data, params, x):
  518. """Call the fit function of the 1D fit model to do the fit.
  519. :param data: The data to fit
  520. :type data: 1D numpy array
  521. :param params: The initial paramters of the fit
  522. :type params: lmfit Parameters
  523. :param x: The data of x axis
  524. :type x: 1D numpy array
  525. :return: The result of the fit
  526. :rtype: lmfit FitResult
  527. """
  528. return self.fitModel.fit(data=data, x=x, params=params, nan_policy='omit')
  529. def _fit_2D(self, data, params, x, y):
  530. """Call the fit function of the 2D fit model to do the fit.
  531. :param data: The flattened data to fit
  532. :type data: 1D numpy array
  533. :param params: The flattened initial paramters of the fit
  534. :type params: lmfit Parameters
  535. :param x: The flattened data of x axis
  536. :type x: 1D numpy array
  537. :param y: The flattened data of y axis
  538. :type y: 1D numpy array
  539. :return: The result of the fit
  540. :rtype: lmfit FitResult
  541. """
  542. data = data.flatten(order='F')
  543. return self.fitModel.fit(data=data, x=x, y=y, params=params, nan_policy='omit')
  544. def fit(self, dataArray, paramsArray, x=None, y=None, input_core_dims=None, dask='parallelized', vectorize=True, keep_attrs=True, daskKwargs=None, **kwargs):
  545. """Call the fit function of the fit model to do the fit.
  546. :param dataArray: The data for the fit
  547. :type dataArray: xarray DataArray
  548. :param paramsArray: The flattened initial paramters of the fit
  549. :type paramsArray: xarray DataArray or lmfit Parameters
  550. :param x: The name of x axis in data or the data of x axis, defaults to None
  551. :type x: str or numpy array, optional
  552. :param y: The name of y axis in data or the data of y axis, defaults to None
  553. :type y: str or numpy array, optional
  554. :param input_core_dims: over write of the same argument in xarray.apply_ufunc, defaults to None, defaults to None
  555. :type input_core_dims: list or array like, optional
  556. :param dask: over write of the same argument in xarray.apply_ufunc, defaults to None, defaults to 'parallelized'
  557. :type dask: str, optional
  558. :param vectorize: over write of the same argument in xarray.apply_ufunc, defaults to None, defaults to True
  559. :type vectorize: bool, optional
  560. :param keep_attrs: over write of the same argument in xarray.apply_ufunc, defaults to None, defaults to True
  561. :type keep_attrs: bool, optional
  562. :param daskKwargs: over write of the same argument in xarray.apply_ufunc, defaults to None, defaults to None
  563. :type daskKwargs: dict, optional
  564. :return: The fit result
  565. :rtype: xarray DataArray
  566. """
  567. kwargs.update(
  568. {
  569. "dask": dask,
  570. "vectorize": vectorize,
  571. "input_core_dims": input_core_dims,
  572. 'keep_attrs': keep_attrs,
  573. }
  574. )
  575. if not daskKwargs is None:
  576. kwargs.update({"dask_gufunc_kwargs": daskKwargs})
  577. if isinstance(paramsArray, type(self.fitModel.make_params())):
  578. if input_core_dims is None:
  579. kwargs.update(
  580. {
  581. "input_core_dims": [['x']],
  582. }
  583. )
  584. if x is None:
  585. if 'x' in dataArray.dims:
  586. x = dataArray['x'].to_numpy()
  587. else:
  588. if isinstance(x, str):
  589. if input_core_dims is None:
  590. kwargs.update(
  591. {
  592. "input_core_dims": [[x]],
  593. }
  594. )
  595. x = dataArray[x].to_numpy()
  596. if self.fitDim == 1:
  597. res = xr.apply_ufunc(self._fit_1D, dataArray, kwargs={'params':paramsArray,'x':x},
  598. output_dtypes=[type(lmfit.model.ModelResult(self.fitModel, self.fitModel.make_params()))],
  599. **kwargs)
  600. if keep_attrs:
  601. res.attrs = copy.copy(dataArray.attrs)
  602. return res
  603. if self.fitDim == 2:
  604. if y is None:
  605. if 'y' in dataArray.dims:
  606. y = dataArray['y'].to_numpy()
  607. if input_core_dims is None:
  608. kwargs.update(
  609. {
  610. "input_core_dims": [['x', 'y']],
  611. }
  612. )
  613. else:
  614. if isinstance(y, str):
  615. kwargs["input_core_dims"][0] = np.append(kwargs["input_core_dims"][0], y)
  616. y = dataArray[y].to_numpy()
  617. elif input_core_dims is None:
  618. kwargs.update(
  619. {
  620. "input_core_dims": [['x', 'y']],
  621. }
  622. )
  623. _x, _y = np.meshgrid(x, y)
  624. _x = _x.flatten()
  625. _y = _y.flatten()
  626. # dataArray = dataArray.stack(_z=(kwargs["input_core_dims"][0][0], kwargs["input_core_dims"][0][1]))
  627. # kwargs["input_core_dims"][0] = ['_z']
  628. res = xr.apply_ufunc(self._fit_2D, dataArray, kwargs={'params':paramsArray,'x':_x, 'y':_y},
  629. output_dtypes=[type(lmfit.model.ModelResult(self.fitModel, self.fitModel.make_params()))],
  630. **kwargs)
  631. if keep_attrs:
  632. res.attrs = copy.copy(dataArray.attrs)
  633. return res
  634. else:
  635. if input_core_dims is None:
  636. kwargs.update(
  637. {
  638. "input_core_dims": [['x'], []],
  639. }
  640. )
  641. if x is None:
  642. if 'x' in dataArray.dims:
  643. x = dataArray['x'].to_numpy()
  644. else:
  645. if isinstance(x, str):
  646. if input_core_dims is None:
  647. kwargs.update(
  648. {
  649. "input_core_dims": [[x], []],
  650. }
  651. )
  652. x = dataArray[x].to_numpy()
  653. if self.fitDim == 1:
  654. return xr.apply_ufunc(self._fit_1D, dataArray, paramsArray, kwargs={'x':x},
  655. output_dtypes=[type(lmfit.model.ModelResult(self.fitModel, self.fitModel.make_params()))],
  656. **kwargs)
  657. if self.fitDim == 2:
  658. if input_core_dims is None:
  659. kwargs.update(
  660. {
  661. "input_core_dims": [['x', 'y'], []],
  662. }
  663. )
  664. if y is None:
  665. if 'y' in dataArray.dims:
  666. y = dataArray['y'].to_numpy()
  667. else:
  668. if isinstance(y, str):
  669. y = dataArray[y].to_numpy()
  670. kwargs["input_core_dims"][0] = np.append(kwargs["input_core_dims"][0], y)
  671. _x, _y = np.meshgrid(x, y)
  672. _x = _x.flatten()
  673. _y = _y.flatten()
  674. # dataArray = dataArray.stack(_z=(kwargs["input_core_dims"][0][0], kwargs["input_core_dims"][0][1]))
  675. # kwargs["input_core_dims"][0] = ['_z']
  676. return xr.apply_ufunc(self._fit_2D, dataArray, paramsArray, kwargs={'x':_x, 'y':_y},
  677. output_dtypes=[type(lmfit.model.ModelResult(self.fitModel, self.fitModel.make_params()))],
  678. **kwargs)
  679. def _eval_1D(self, fitResult, x):
  680. """Call the eval function of the 1D fit model to calculate the curve.
  681. :param fitResult: The result of fit
  682. :type fitResult: lmfit FitResult
  683. :param x: The data of x axis
  684. :type x: 1D numpy array
  685. :return: The curve based on the fit result
  686. :rtype: 1D numpy array
  687. """
  688. return self.fitModel.eval(x=x, **fitResult.best_values)
  689. def _eval_2D(self, fitResult, x, y, shape):
  690. """Call the eval function of the 2D fit model to calculate the curve.
  691. :param fitResult: The result of fit
  692. :type fitResult: lmfit FitResult
  693. :param x: The data of x axis
  694. :type x: 1D numpy array
  695. :param y: The flattened data of y axis
  696. :type y: 1D numpy array
  697. :param shape: The desired shape for output
  698. :type shape: list, tuple or array like
  699. :return: The curve based on the fit result
  700. :rtype: 2D numpy array
  701. """
  702. res = self.fitModel.eval(x=x, y=y, **fitResult.best_values)
  703. return res.reshape(shape, order='F')
  704. def eval(self, fitResultArray, x=None, y=None, output_core_dims=None, prefix="", dask='parallelized', vectorize=True, daskKwargs=None, **kwargs):
  705. """Call the eval function of the fit model to calculate the curve.
  706. :param fitResultArray: The result of fit
  707. :type fitResultArray: xarray DataArray
  708. :param x: The name of x axis in data or the data of x axis, defaults to None
  709. :type x: str or numpy array, optional
  710. :param y: The name of y axis in data or the data of y axis, defaults to None
  711. :type y: str or numpy array, optional
  712. :param output_core_dims: over write of the same argument in xarray.apply_ufunc, defaults to None
  713. :type output_core_dims: list, optional
  714. :param prefix: prefix str of the fit model, defaults to ""
  715. :type prefix: str, optional
  716. :param dask: over write of the same argument in xarray.apply_ufunc, defaults to 'parallelized'
  717. :type dask: str, optional
  718. :param vectorize: over write of the same argument in xarray.apply_ufunc, defaults to True
  719. :type vectorize: bool, optional
  720. :param daskKwargs: over write of the same argument in xarray.apply_ufunc, defaults to None
  721. :type daskKwargs: dict, optional
  722. :return: The curve based on the fit result
  723. :rtype: xarray
  724. """
  725. kwargs.update(
  726. {
  727. "dask": dask,
  728. "vectorize": vectorize,
  729. "output_core_dims": output_core_dims,
  730. }
  731. )
  732. if daskKwargs is None:
  733. daskKwargs = {}
  734. if self.fitDim == 1:
  735. if output_core_dims is None:
  736. kwargs.update(
  737. {
  738. "output_core_dims": prefix+'x',
  739. }
  740. )
  741. output_core_dims = [prefix+'x']
  742. daskKwargs.update(
  743. {
  744. 'output_sizes': {
  745. output_core_dims[0]: np.size(x),
  746. },
  747. 'meta': np.ndarray((0,0), dtype=float)
  748. }
  749. )
  750. kwargs.update(
  751. {
  752. "dask_gufunc_kwargs": daskKwargs,
  753. }
  754. )
  755. res = xr.apply_ufunc(self._eval_1D, fitResultArray, kwargs={"x":x}, **kwargs)
  756. return res.assign_coords({prefix+'x':np.array(x)})
  757. if self.fitDim == 2:
  758. if output_core_dims is None:
  759. kwargs.update(
  760. {
  761. "output_core_dims": [[prefix+'x', prefix+'y']],
  762. }
  763. )
  764. output_core_dims = [prefix+'x', prefix+'y']
  765. daskKwargs.update(
  766. {
  767. 'output_sizes': {
  768. output_core_dims[0]: np.size(x),
  769. output_core_dims[1]: np.size(y),
  770. },
  771. 'meta': np.ndarray((0,0), dtype=float)
  772. },
  773. )
  774. kwargs.update(
  775. {
  776. "dask_gufunc_kwargs": daskKwargs,
  777. }
  778. )
  779. _x, _y = np.meshgrid(x, y)
  780. _x = _x.flatten()
  781. _y = _y.flatten()
  782. res = xr.apply_ufunc(self._eval_2D, fitResultArray, kwargs={"x":_x, "y":_y, "shape":(len(x), len(y))}, **kwargs)
  783. return res.assign_coords({prefix+'x':np.array(x), prefix+'y':np.array(y)})
  784. def _get_fit_value_single(self, fitResult, key):
  785. """get value of one single parameter from fit result
  786. :param fitResult: The result of the fit
  787. :type fitResult: lmfit FitResult
  788. :param key: The name of the parameter
  789. :type key: str
  790. :return: The vaule of the parameter
  791. :rtype: float
  792. """
  793. return fitResult.params[key].value
  794. def _get_fit_value(self, fitResult, params):
  795. """get value of parameters from fit result
  796. :param fitResult: The result of the fit
  797. :type fitResult: lmfit FitResult
  798. :param params: A list of the name of paramerters to read
  799. :type params: list or array like
  800. :return: The vaule of the parameter
  801. :rtype: 1D numpy array
  802. """
  803. func = np.vectorize(self._get_fit_value_single)
  804. res = tuple(
  805. func(fitResult, key)
  806. for key in params
  807. )
  808. return res
  809. def get_fit_value(self, fitResult, dask='parallelized', **kwargs):
  810. """get value of parameters from fit result
  811. :param fitResult: The result of the fit
  812. :type fitResult: lmfit FitResult
  813. :param dask: over write of the same argument in xarray.apply_ufunc, defaults to 'parallelized'
  814. :type dask: str, optional
  815. :return: The vaule of the parameter
  816. :rtype: xarray DataArray
  817. """
  818. firstIndex = {
  819. key: fitResult[key][0]
  820. for key in fitResult.dims
  821. }
  822. firstFitResult = fitResult.sel(firstIndex).item()
  823. params = list(firstFitResult.params.keys())
  824. output_core_dims=[ [] for _ in range(len(params))]
  825. kwargs.update(
  826. {
  827. "dask": dask,
  828. "output_core_dims": output_core_dims,
  829. }
  830. )
  831. value = xr.apply_ufunc(self._get_fit_value, fitResult, kwargs=dict(params=params), **kwargs)
  832. value = xr.Dataset(
  833. data_vars={
  834. params[i]: value[i]
  835. for i in range(len(params))
  836. },
  837. attrs=fitResult.attrs
  838. )
  839. return value
  840. def _get_fit_std_single(self, fitResult, key):
  841. """get standard deviation of one single parameter from fit result
  842. :param fitResult: The result of the fit
  843. :type fitResult: lmfit FitResult
  844. :param key: The name of the parameter
  845. :type key: str
  846. :return: The vaule of the parameter
  847. :rtype: float
  848. """
  849. return fitResult.params[key].stderr
  850. def _get_fit_std(self, fitResult, params):
  851. """get standard deviation of parameters from fit result
  852. :param fitResult: The result of the fit
  853. :type fitResult: lmfit FitResult
  854. :param params: A list of the name of paramerters to read
  855. :type params: list or array like
  856. :return: The vaule of the parameter
  857. :rtype: 1D numpy array
  858. """
  859. func = np.vectorize(self._get_fit_std_single)
  860. res = tuple(
  861. func(fitResult, key)
  862. for key in params
  863. )
  864. return res
  865. def get_fit_std(self, fitResult, dask='parallelized', **kwargs):
  866. """get standard deviation of parameters from fit result
  867. :param fitResult: The result of the fit
  868. :type fitResult: lmfit FitResult
  869. :param dask: over write of the same argument in xarray.apply_ufunc, defaults to 'parallelized'
  870. :type dask: str, optional
  871. :return: The vaule of the parameter
  872. :rtype: xarray DataArray
  873. """
  874. firstIndex = {
  875. key: fitResult[key][0]
  876. for key in fitResult.dims
  877. }
  878. firstFitResult = fitResult.sel(firstIndex).item()
  879. params = list(firstFitResult.params.keys())
  880. output_core_dims=[ [] for _ in range(len(params))]
  881. kwargs.update(
  882. {
  883. "dask": dask,
  884. "output_core_dims": output_core_dims,
  885. }
  886. )
  887. value = xr.apply_ufunc(self._get_fit_std, fitResult, kwargs=dict(params=params), **kwargs)
  888. value = xr.Dataset(
  889. data_vars={
  890. params[i]: value[i]
  891. for i in range(len(params))
  892. },
  893. attrs=fitResult.attrs
  894. )
  895. return value
  896. def _get_fit_full_result_single(self, fitResult, key):
  897. """get the full result of one single parameter from fit result
  898. :param fitResult: The result of the fit
  899. :type fitResult: lmfit FitResult
  900. :param key: The name of the parameter
  901. :type key: str
  902. :return: The vaule of the parameter
  903. :rtype: float
  904. """
  905. if not fitResult.params[key].value is None:
  906. value = fitResult.params[key].value
  907. else:
  908. value = np.nan
  909. if not fitResult.params[key].stderr is None:
  910. std = fitResult.params[key].stderr
  911. else:
  912. std = np.nan
  913. return ufloat(value, std)
  914. def _get_fit_full_result(self, fitResult, params):
  915. """get the full result of parameters from fit result
  916. :param fitResult: The result of the fit
  917. :type fitResult: lmfit FitResult
  918. :param params: A list of the name of paramerters to read
  919. :type params: list or array like
  920. :return: The vaule of the parameter
  921. :rtype: 1D numpy array
  922. """
  923. func = np.vectorize(self._get_fit_full_result_single)
  924. res = tuple(
  925. func(fitResult, key)
  926. for key in params
  927. )
  928. return res
  929. def get_fit_full_result(self, fitResult, dask='parallelized', **kwargs):
  930. """get the full result of parameters from fit result
  931. :param fitResult: The result of the fit
  932. :type fitResult: lmfit FitResult
  933. :param dask: over write of the same argument in xarray.apply_ufunc, defaults to 'parallelized'
  934. :type dask: str, optional
  935. :return: The vaule of the parameter
  936. :rtype: xarray DataArray
  937. """
  938. firstIndex = {
  939. key: fitResult[key][0]
  940. for key in fitResult.dims
  941. }
  942. firstFitResult = fitResult.sel(firstIndex).item()
  943. params = list(firstFitResult.params.keys())
  944. output_core_dims=[ [] for _ in range(len(params))]
  945. kwargs.update(
  946. {
  947. "dask": dask,
  948. "output_core_dims": output_core_dims,
  949. }
  950. )
  951. value = xr.apply_ufunc(self._get_fit_full_result, fitResult, kwargs=dict(params=params), **kwargs)
  952. value = xr.Dataset(
  953. data_vars={
  954. params[i]: value[i]
  955. for i in range(len(params))
  956. },
  957. attrs=fitResult.attrs
  958. )
  959. return value