calculation of moessbauer spectra More...
#include <Moesscalc.h>
| Public Member Functions | |
| void | calc () | 
| fill Hamiltonian, solve and precalculate transition propabilities | |
| double | get_Egamma (int f, int i) | 
| get Transition Energy from initial state i to final State f | |
| double | Vzz_mms (double Vzz) | 
| transform V/Angstroem^2 into mm/s at axialsymmetric field gradient | |
| double | get_dv_Douplett () | 
| get current quadrupol splitting in mm/s | |
| double | I (double beta, double gamma, int f, int i) | 
| Transition propability from i to f in crystal geometry.  More... | |
| void | I (double beta, double gamma, double alpha, int f, int i, double *Int, double *phi, double *xi) | 
| Transition propability of polarized radiation in crystal geometry.  More... | |
| double | Idiag (double beta, int f, int i) | 
| transition propabilities for a diagonal hamitonian  More... | |
| double | I (int f, int i) | 
| powder transition propabilities for the current hamiltonian | |
| double | f_lorenz (double x, double dx, double x0, double I) | 
| Lorentz functions I/(pi*dx)*1/(1+((x-x0)/dx)^2) | |
| double | f_lorenz_fast (double x, double invdx, double x0, double I) | 
| Lorentz functions I*invdx/((1+((x-x0)*invdx)^2)*M_PI) | |
| void | get_SC_spectrum (double *I, int cc, double *v, double *p) | 
| single crystal spectrum, static Hamiltonian | |
| void | get_Powder_spectrum (double *I, int cc, double *v, double *p) | 
| powder spectrum, static Hamiltonian | |
| void | get_SrcFld_spectrum (double *I, int cc, double *v, double *p) | 
| static Hamiltonian with source in field | |
| void | get_fm_Powder_spectrum (double *I, int cc, double *v, double *p) | 
| ferromagnet in transverse field spectrum | |
| void | get_afm_Powder_spectrum (double *I, int cc, double *v, double *p) | 
| antiferromagnet in transverse field spetrum | |
| void | get_PseudoVoigt_profile (double *I, int cc, double *v, double *p) | 
| void | get_FeCal_spectrum (double *I, int cc, double *v, double *p, double U) | 
| iron calibrations, recalculates sin-velocity using monitur voltage U | |
| void | get_Dyn_spectrum (double *I, int cc, double *v, double *p) | 
| Blume-Jones model with arbitrary number of intermediate states.  More... | |
| void | get_MLR_spectrum (double *I, int cc, double *v, double *p) | 
| multi level relaxation model, asuming a potential barrier  More... | |
| void | get_DynDiagDiag_spectrum (double *I, int cc, double *v, double *p) | 
| Fast Blume Jones for two states with diagonal Hamiltonian. | |
| void | get_DynDiagDiagDiag_spectrum (double *I, int cc, double *v, double *p) | 
| Fast Blume Jones for three states with diagonal Hamiltonian. | |
| void | add_noise (double *I, int cc) | 
| adds gaussian noise to spectrum I with standard deviation sqrt(I) | |
| void | get_FeCaltriang_spectrum (double *I, int cc, double *v, double *p, double U) | 
| iron calibrations, recalculates triangular-velocity | |
| void | get_FeAs_spectrum (double *I, int cc, double *v, double *p) | 
| calculates FeAs spectrum (Europhys. Lett., 9 (l), pp. 87-92 (1989)) | |
| void | SimpleInvert (double *x, double *y, int count, int steps, double **Mspec, double *rho) | 
| Hesse Ruebartsch without smoothing parameter.  More... | |
| 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  More... | |
| 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.  More... | |
| void | HesseRuebartsch (double **M0, double *rho, double *y, int steps, int cc, double smooth) | 
| Hesse Ruebartsch parameter distribution.  More... | |
| void | LeCear (double **M0, double *rho, double *y, int steps, int cc, double smooth) | 
| Hesse Ruebartsch, but restricted to positive weights.  More... | |
| double | dv_Debye (double T, double TD, double Meff) | 
| quadratic dopplereffect in Debye approximation | |
| double | Abs_Debye (double T, double TD, double Meff) | 
| Debye Waller factor in Debye approximation. | |
| double | sWave_SuperfluidDensity (double SF0, double Tc, double T, double delta0) | 
| superfluid density in s-wave symmetry  More... | |
| double | Hc2_WHH (double T, double Tc, double dHc2dT, double alpha, double lambdaSO) | 
| Werthamer-Helfand-Honenberg (WHH) model for upper critial field.  More... | |
| double | OrderParameterDistr (double T, double Tn, double sigma, double alpha, double beta) | 
| order parameter temperature dependence  More... | |
| double | errorfunction (double arg) | 
| Gaussian error function f(x)=2/sqrt(pi)*int(0,x)[exp(-a^2)]da. | |
| 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.  More... | |
| void | afm_flipflop (double Bex, double BJ, double BA, double zeta, double *theta1, double *theta2) | 
| calculates tilting of antiferromagnetic (afm) moments in field  More... | |
| void | fm_flip (double Bex, double BA, double zeta, double *theta1) | 
| see afm_fliflop, but without antiferromagnetic coupling | |
calculation of moessbauer spectra
spectrum constructors of the form spectrum(I,cc,v,p) need the following arguments
| I,: | intensity, array of length cc for output | 
| cc,: | channel count | 
| v,: | velocity array of length cc | 
| parameters array | 
| void MoessCalc::afm_flipflop | ( | double | Bex, | 
| double | BJ, | ||
| double | BA, | ||
| double | zeta, | ||
| double * | theta1, | ||
| double * | theta2 | ||
| ) | 
calculates tilting of antiferromagnetic (afm) moments in field
| Bex,: | external field | 
| BJ,: | afm coupling constant | 
| BA,: | anisotropy field | 
| zeta,: | tilting of anisotropy axis with respect to Bex | 
| theta1,: | polar angle of moment 1 to Bex | 
| theta2,: | polar angle of moment 2 to Bex | 
| double MoessCalc::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
simulates the spectrum with a high number of Lorentzian singlets using SimpleInvert function
| void MoessCalc::get_Dyn_spectrum | ( | double * | I, | 
| int | cc, | ||
| double * | v, | ||
| double * | p | ||
| ) | 
Blume-Jones model with arbitrary number of intermediate states.
Physical Review 174.2 (1968): 351
| void MoessCalc::get_MLR_spectrum | ( | double * | I, | 
| int | cc, | ||
| double * | v, | ||
| double * | p | ||
| ) | 
multi level relaxation model, asuming a potential barrier
Journal of Physics: Condensed Matter 23 (2011), Nr. 42, S.426003
| double MoessCalc::Hc2_WHH | ( | double | T, | 
| double | Tc, | ||
| double | dHc2dT, | ||
| double | alpha, | ||
| double | lambdaSO | ||
| ) | 
Werthamer-Helfand-Honenberg (WHH) model for upper critial field.
see Collings, EW: Applied Superconductivity, Metallurgy, and Physics of Titanium Alloys: Volume 1: Fundamentals. Bd. 1. Springer 1986
| void MoessCalc::HesseRuebartsch | ( | double ** | M0, | 
| double * | rho, | ||
| double * | y, | ||
| int | steps, | ||
| int | cc, | ||
| double | smooth | ||
| ) | 
Hesse Ruebartsch parameter distribution.
(1974 J. Phys. E: Sci. Instrum. 7 526), uses Linear Algebra to weihght a number of subspectra to meet both best fitting and a smooth weight histogram
| M0,: | Matrix of subspectra | 
| rho,: | spectra weights to be find | 
| y,: | data points | 
| steps,: | subspectra count | 
| cc,: | channel count of each subspectrum | 
| smooth,: | smoothing parameter | 
| double MoessCalc::I | ( | double | beta, | 
| double | gamma, | ||
| int | f, | ||
| int | i | ||
| ) | 
Transition propability from i to f in crystal geometry.
| beta,: | polar angle of gamma beam to field gradient z-axis | 
| gamma,: | azimutal angle of gamma beam to field gradient z-axis | 
| f,: | final state index | 
| i,: | initial state index | 
| void MoessCalc::I | ( | double | beta, | 
| double | gamma, | ||
| double | alpha, | ||
| int | f, | ||
| int | i, | ||
| double * | Int, | ||
| double * | phi, | ||
| double * | xi | ||
| ) | 
Transition propability of polarized radiation in crystal geometry.
| beta,: | polar angle of gamma beam to field gradient z-axis | 
| gamma,: | azimutal angle of gamma beam to field gradient z-axis | 
| alpha,: | azimutal angle between source and sample z-axis | 
| f,: | final state index | 
| i,: | initial state index | 
| Int,: | intensity output | 
| phi,: | tilting of the polarisation ellipse with respect to z-axis | 
| xi,: | eccentricity of the polarization matrix | 
| double MoessCalc::Idiag | ( | double | beta, | 
| int | f, | ||
| int | i | ||
| ) | 
transition propabilities for a diagonal hamitonian
| beta,: | incoming polar angle | 
| f,: | finale state index | 
| i,: | initial state index | 
| void MoessCalc::LeCear | ( | double ** | M0, | 
| double * | rho, | ||
| double * | y, | ||
| int | steps, | ||
| int | cc, | ||
| double | smooth | ||
| ) | 
Hesse Ruebartsch, but restricted to positive weights.
Journal of Physics E: Scientifc Instruments 12 (1979), Nr. 11, S. 1083 
Instead of using linear algebra the determination of the weight histogram might be achieved by numeric minimization. For parameters see HesseRuebartsch 
| void MoessCalc::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.
ambridge algorithm: Mon. Not. R Astr. Soc. (1984) 211, 111-124
 MEM in Moessbauer: Brand, R. A., and G. Le Caër. "Improving the validity of Mössbauer hyperfine parameter distributions: the maximum entropy formalism and its applications."
| M0,: | subspectra to be weighted | 
| rho,: | weight histogram to be smoothed | 
| rhoerr,: | errors of the weights | 
| y,: | experimental spectrum | 
| I0,: | base line | 
| steps,: | subspectra counts | 
| cc,: | channel count | 
| stepsize,: | of algorithm | 
| MaxLambda,: | maximum fitting weight | 
| MaxDif,: | maximum accepted divergence off entropy gradient and chi2 gradient | 
| double MoessCalc::OrderParameterDistr | ( | double | T, | 
| double | Tn, | ||
| double | sigma, | ||
| double | alpha, | ||
| double | beta | ||
| ) | 
order parameter temperature dependence
with two critical exponents and a gaussian distributed transition temperature:
 OP(T)=int(0,)[Gauss(t,Tn,sigma)*(1-(T/t)^alpha)^beta]dt
| T,: | temperature | 
| Tn,: | transition temperature mean | 
| sigma,: | standarrd deviation of the transition temperature | 
| alpha,: | first critical exponent | 
| beta,: | second critical exponent | 
| void MoessCalc::SimpleInvert | ( | double * | x, | 
| double * | y, | ||
| int | count, | ||
| int | steps, | ||
| double ** | Mspec, | ||
| double * | rho | ||
| ) | 
Hesse Ruebartsch without smoothing parameter.
a number of subspectra are weighted to fit the total spectrum in the best way using linear algebra in the same way as Hesse-Ruebartsch, but without smoothing
| x,: | velocity data | 
| y,: | intensity data | 
| count,: | channel count of each spectrum | 
| steps,: | subspectra count | 
| Mspec,: | count*steps Matrix with subspectra | 
| rho,: | weights of the subspectra | 
| double MoessCalc::sWave_SuperfluidDensity | ( | double | SF0, | 
| double | Tc, | ||
| double | T, | ||
| double | delta0 | ||
| ) | 
superfluid density in s-wave symmetry
Supercond. Sci. Technol.21(2008) 082003 (3pp)
| SF0,: | superfluid density at T=0 | 
| Tc,: | critical temperature | 
| T,: | temperature | 
| delta0,: | gap value | 
| void MoessCalc::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.
acording to the transmission integral the final absorption spectrum is the convolution of the Lorentzian source line with the absorber spectrum. The absorber spectrum represents the crosssection for i.e. it appears in the exponent of an exponential function.
 A symmetric velocity specing considerably can reduce the calculation costs for the convolution process. However the function can deal with asymmetric spacing as well.
| SHI,: | thin absorber spectrum | 
| v,: | velocities | 
| count,: | channel count | 
| A,: | spectral area of the thin-absorber-spectrum | 
| I0,: | baseline | 
| ta,: | effective thickness | 
| omegasrc,: | source line width | 
| omegaabs,: | absorber line width | 
| fr,: | resonant fraction, i.e. the part of the photon count rate which in principle is able to be resonantly absorbed, i.e. it's not background radiation, it hits the sample, and it respects both source and sample Debye Waller factor | 
 1.8.5
 1.8.5