My Project
 All Classes Functions Variables Pages
Moesscalc.h
1  #ifndef MOESSCALC_H
2 #define MOESSCALC_H
3 #define DegRad 0.01745329251994329576923690768489
4 #define RadDeg 57.295779513082320876798154814105
5 
7 {
8 double a, b, c;
9 };
10 
11 /* Calculation of energy levels and transition propabilities:
12  all values in mm/s!
13 */
14 #include <QtGlobal>
15 #include <QMessageBox>
16 #include <QProgressDialog>
17 #include <gsl/gsl_vector.h>
18 #include <gsl/gsl_complex.h>
19 #include <gsl/gsl_complex_math.h>
20 #include <gsl/gsl_eigen.h>
21 #include <gsl/gsl_integration.h>
22 #include <gsl/gsl_math.h>
23 #include <gsl/gsl_matrix.h>
24 #include <gsl/gsl_linalg.h>
25 #include <gsl/gsl_blas.h>
26 #include <gsl/gsl_multifit_nlin.h>
27 #include <gsl/gsl_multimin.h>
28 #include <gsl/gsl_multiroots.h>
29 #include <gsl/gsl_roots.h>
30 #include <gsl/gsl_sf.h>
31 #include <gsl/gsl_sort_vector.h>
32 #include <gsl/gsl_cdf.h>
33 #include <gsl/gsl_rng.h>
34 #include <gsl/gsl_const_mksa.h>
35 #include <stdlib.h>
36 #include "plot2D.h"
37 #include <cmath>
38 #include <QTime>
39 
40 
42 typedef struct LeCearContainer{
43  gsl_matrix *M;
44  gsl_matrix *D0;
45  gsl_matrix *D1;
46  gsl_matrix *A;
47  gsl_vector *S;
48  int steps;
49  double lambda;
50  double count;
52 
54 enum alphamode{decrease_in,decrease_out,increase_out, increase_in,
55  initialize, exceed_aim};
56 
58 
65 class MoessCalc{
66  public:
67  MoessCalc();
68  ~MoessCalc();
70  void calc();
72  double get_Egamma(int f, int i);
74  double Vzz_mms(double Vzz);
76  double get_dv_Douplett();
78 
83  double I(double beta, double gamma, int f, int i);
85 
94  void I(double beta, double gamma, double alpha, int f, int i,
95  double* Int, double* phi, double* xi);
97 
101  double Idiag(double beta, int f, int i);
103  double I(int f, int i);
105  double f_lorenz(double x, double dx,double x0,double I);
107  double f_lorenz_fast(double x, double invdx,double x0,double I);
109  void get_SC_spectrum(double *I, int cc , double* v, double* p);
111  void get_Powder_spectrum(double *I, int cc , double* v, double* p);
113  void get_SrcFld_spectrum(double *I, int cc , double* v, double* p);
115  void get_fm_Powder_spectrum(double *I, int cc , double* v, double* p);
117  void get_afm_Powder_spectrum(double *I, int cc , double* v, double* p);
118  void get_PseudoVoigt_profile(double *I, int cc , double* v, double* p);
120  void get_FeCal_spectrum(double *I, int cc , double* v, double* p,
121  double U);
123 
124  void get_Dyn_spectrum(double *I, int cc , double* v, double* p);
126 
127  void get_MLR_spectrum(double *I, int cc , double* v, double* p);
129  void get_DynDiagDiag_spectrum(double *I, int cc , double* v, double* p);
131  void get_DynDiagDiagDiag_spectrum(double *I, int cc , double* v, double* p);
133  void add_noise(double *I, int cc);
135  void get_FeCaltriang_spectrum(double *I, int cc , double* v, double* p,
136  double U);
138  void get_FeAs_spectrum(double *I, int cc , double* v, double* p);
140 
150  void SimpleInvert(double* x, double *y, int count, int steps,
151  double **Mspec, double* rho);
153 
156  double baseline(double* x,double* y,int count,
157  double* Baseline,double* Face,
158  double bperx, double omega,
159  double** rho,double** rhox, int *rhosteps, double* INew);
161 
175  void MEM_Cambridge(double**M0, double* rho, double* rhoerr, double *y,
176  double I0, int steps, int cc, double stepsize,
177  int MaxIterations, double MaxLambda,
178  double MaxDif, double status[5]);
180 
190  void HesseRuebartsch(double**M0, double* rho, double *y,
191  int steps, int cc, double smooth);
193 
198  void LeCear(double**M0, double* rho, double *y,
199  int steps, int cc, double smooth);
201  double dv_Debye(double T, double TD, double Meff);
203  double Abs_Debye(double T, double TD, double Meff);
205 
211  double sWave_SuperfluidDensity(double SF0, double Tc, double T,
212  double delta0);
214 
218  double Hc2_WHH(double T, double Tc, double dHc2dT, double alpha,
219  double lambdaSO);
221 
230  double OrderParameterDistr(double T, double Tn, double sigma,
231  double alpha, double beta);
233  double errorfunction(double arg);
235 
255  void TransmissionIntegral(double *SHI,double *v, int count,
256  double A, double I0, double ta,
257  double omegasrc, double omegaabs, double fr);
259 
266  void afm_flipflop(double Bex, double BJ, double BA, double zeta,
267  double* theta1, double*theta2);
269  void fm_flip(double Bex, double BA, double zeta, double* theta1);
270  private:
272  void fill_V();
273  void set_parameters(double B, double Vzz,double theta, double phi,
274  double eta);
275  gsl_matrix_complex *V0;
276  gsl_matrix_complex *V1;
277  gsl_matrix_complex *C0;
278  gsl_matrix_complex *C1;
279 
282  void SolveTridiagSys(gsl_matrix_complex* M, gsl_vector_complex* X,
283  gsl_vector_complex* b);
284 
285  double B;
286  double Vzz;
287  double theta;
288  double phi;
289  double eta;
290  double I0;
291  double E0[2],E1[4];
292 
294  double a;
296 
297  double b;
298 
300  double b0;
301  int firstcalc;
302 
303  void init();
304 
307  double c,mun,Egamma,Q,g12,g32,Runi,Erec,kbolz,hbar_c_over_Egamma;
309 
311 
314  void T_sigma(double *SHI,double *v, int count,
315  double A, double I0, double* OneMinusT,
316  double ta, double omegasrc, double omegaabs);
317 
319 
322  int steps;
324  int count;
325  gsl_vector *dchi2;
326  gsl_matrix *d2chi2;
327  gsl_vector *dS;
328  gsl_matrix *d2S;
329  gsl_vector *yS;
330  gsl_vector *chi;
331  gsl_vector *P;
332  gsl_vector *vtemp;
333  gsl_matrix *M;
334  double cangle;
335  double cdif;
336  double lagrange;
337  double Pstepmax;
338  void MEM_init(double**M, double *y);
339  void MEM_free();
340  void MEM_calc_derivatives();
341  void MEM_calc_diff();
342  // now specific Cambridge algorithm content follows
343  int ssEWc;
344  gsl_vector *e[4];
345  gsl_vector *eON[4];
346  gsl_vector *ssEW;
347  gsl_matrix *g;
348  gsl_matrix *ssg;
349  gsl_matrix *ssEV;
350  gsl_matrix *ssM;
351  gsl_matrix *MEV;
352  gsl_vector *MEW;
353  void MEM_Cambridge_create_subspace();
354  void MEM_Cambridge_free_subspace();
355  void MEM_Cambridge_calc_step();
356  void MEM_error_analysis(double* rho, double *rhoerror, double lambda);
357  double MEM_Cambridge_lamda();
358  double MEM_Cambridge_chi2();
359  void MEM_Cambridge_tune_lambda(double* rho, double MaxLambda,
361  double status[5]);
363 
364 
366 
373  void fitParaboloid(double** r, double* y, int N,int dim,
374  double* result, double* errors);
376  double sum(double **x, int N, int dim, int* indeces);
377 
378 
381  gsl_eigen_hermv_workspace * ws1;
382  gsl_vector *eval1;
383  gsl_eigen_hermv_workspace * ws0;
384  gsl_vector *eval0;
386 
389  gsl_complex c12m12[2];
390  gsl_complex c12p12[2];
391  gsl_complex c32m32[4];
392  gsl_complex c32m12[4];
393  gsl_complex c32p12[4];
394  gsl_complex c32p32[4];
396 
399  gsl_complex Km[2][4];
400  gsl_complex K0[2][4];
401  gsl_complex Kp[2][4];
403 };
406 double Hafm_f(const gsl_vector * x, void * params );
407 void Hafm_df(const gsl_vector * x, void * params, gsl_vector * g );
408 void Hafm_fdf(const gsl_vector * x, void * params,double * f, gsl_vector *g);
410 
413 double Hfm_f(const gsl_vector * x, void * params );
414 void Hfm_df(const gsl_vector * x, void * params, gsl_vector * g );
415 void Hfm_fdf(const gsl_vector * x, void * params,double * f, gsl_vector *g);
417 
418 #endif // MOESSCALC_H
double f_lorenz(double x, double dx, double x0, double I)
Lorentz functions I/(pi*dx)*1/(1+((x-x0)/dx)^2)
Definition: Moesscalc.cpp:123
void LeCear(double **M0, double *rho, double *y, int steps, int cc, double smooth)
Hesse Ruebartsch, but restricted to positive weights.
Definition: Moesscalc.cpp:371
void calc()
fill Hamiltonian, solve and precalculate transition propabilities
Definition: Moesscalc.cpp:36
Definition: Moesscalc.h:6
double errorfunction(double arg)
Gaussian error function f(x)=2/sqrt(pi)*int(0,x)[exp(-a^2)]da.
Definition: Moesscalc.cpp:2512
void MEM_Cambridge(double **M0, double *rho, double *rhoerr, double *y, double I0, int steps, int cc, double stepsize, int MaxIterations, double MaxLambda, double MaxDif, double status[5])
Maximum Entroy method using Cambridge algorithm to track MEM path.
Definition: Moesscalc.cpp:761
double Vzz_mms(double Vzz)
transform V/Angstroem^2 into mm/s at axialsymmetric field gradient
Definition: Moesscalc.cpp:1125
calculation of moessbauer spectra
Definition: Moesscalc.h:65
void get_FeAs_spectrum(double *I, int cc, double *v, double *p)
calculates FeAs spectrum (Europhys. Lett., 9 (l), pp. 87-92 (1989))
Definition: Moesscalc.cpp:2557
double Abs_Debye(double T, double TD, double Meff)
Debye Waller factor in Debye approximation.
Definition: Moesscalc.cpp:2334
void get_FeCal_spectrum(double *I, int cc, double *v, double *p, double U)
iron calibrations, recalculates sin-velocity using monitur voltage U
Definition: Moesscalc.cpp:1808
double get_dv_Douplett()
get current quadrupol splitting in mm/s
Definition: Moesscalc.cpp:167
void get_afm_Powder_spectrum(double *I, int cc, double *v, double *p)
antiferromagnet in transverse field spetrum
Definition: Moesscalc.cpp:1554
void get_fm_Powder_spectrum(double *I, int cc, double *v, double *p)
ferromagnet in transverse field spectrum
Definition: Moesscalc.cpp:1479
void fm_flip(double Bex, double BA, double zeta, double *theta1)
see afm_fliflop, but without antiferromagnetic coupling
Definition: Moesscalc.cpp:1636
void get_SC_spectrum(double *I, int cc, double *v, double *p)
single crystal spectrum, static Hamiltonian
Definition: Moesscalc.cpp:1386
void afm_flipflop(double Bex, double BJ, double BA, double zeta, double *theta1, double *theta2)
calculates tilting of antiferromagnetic (afm) moments in field
Definition: Moesscalc.cpp:1669
void get_FeCaltriang_spectrum(double *I, int cc, double *v, double *p, double U)
iron calibrations, recalculates triangular-velocity
Definition: Moesscalc.cpp:1843
void add_noise(double *I, int cc)
adds gaussian noise to spectrum I with standard deviation sqrt(I)
Definition: Moesscalc.cpp:320
void get_DynDiagDiagDiag_spectrum(double *I, int cc, double *v, double *p)
Fast Blume Jones for three states with diagonal Hamiltonian.
Definition: Moesscalc.cpp:2233
void get_MLR_spectrum(double *I, int cc, double *v, double *p)
multi level relaxation model, asuming a potential barrier
Definition: Moesscalc.cpp:2102
double f_lorenz_fast(double x, double invdx, double x0, double I)
Lorentz functions I*invdx/((1+((x-x0)*invdx)^2)*M_PI)
Definition: Moesscalc.cpp:127
void get_Dyn_spectrum(double *I, int cc, double *v, double *p)
Blume-Jones model with arbitrary number of intermediate states.
Definition: Moesscalc.cpp:1878
double Idiag(double beta, int f, int i)
transition propabilities for a diagonal hamitonian
Definition: Moesscalc.cpp:306
double sWave_SuperfluidDensity(double SF0, double Tc, double T, double delta0)
superfluid density in s-wave symmetry
Definition: Moesscalc.cpp:2357
double OrderParameterDistr(double T, double Tn, double sigma, double alpha, double beta)
order parameter temperature dependence
Definition: Moesscalc.cpp:2526
void get_DynDiagDiag_spectrum(double *I, int cc, double *v, double *p)
Fast Blume Jones for two states with diagonal Hamiltonian.
Definition: Moesscalc.cpp:2023
double I(double beta, double gamma, int f, int i)
Transition propability from i to f in crystal geometry.
Definition: Moesscalc.cpp:170
void get_SrcFld_spectrum(double *I, int cc, double *v, double *p)
static Hamiltonian with source in field
Definition: Moesscalc.cpp:1416
double get_Egamma(int f, int i)
get Transition Energy from initial state i to final State f
Definition: Moesscalc.cpp:300
double baseline(double *x, double *y, int count, double *Baseline, double *Face, double bperx, double omega, double **rho, double **rhox, int *rhosteps, double *INew)
determination of baseline and absorption area
Definition: Moesscalc.cpp:571
void HesseRuebartsch(double **M0, double *rho, double *y, int steps, int cc, double smooth)
Hesse Ruebartsch parameter distribution.
Definition: Moesscalc.cpp:478
void TransmissionIntegral(double *SHI, double *v, int count, double A, double I0, double ta, double omegasrc, double omegaabs, double fr)
Convolution of source line with absorber spectrum.
Definition: Moesscalc.cpp:2381
void SimpleInvert(double *x, double *y, int count, int steps, double **Mspec, double *rho)
Hesse Ruebartsch without smoothing parameter.
Definition: Moesscalc.cpp:1329
void get_Powder_spectrum(double *I, int cc, double *v, double *p)
powder spectrum, static Hamiltonian
Definition: Moesscalc.cpp:1779
parameters to execute GSL fitting of parameter distribution
Definition: Moesscalc.h:42
double Hc2_WHH(double T, double Tc, double dHc2dT, double alpha, double lambdaSO)
Werthamer-Helfand-Honenberg (WHH) model for upper critial field.
Definition: Moesscalc.cpp:2470
double dv_Debye(double T, double TD, double Meff)
quadratic dopplereffect in Debye approximation
Definition: Moesscalc.cpp:2328