DMRITool  v0.1.1-139-g860d86b4
Diffusion MRI Tool
itkSpecialFunctionGenerator.h
Go to the documentation of this file.
1 
18 #ifndef __itkSpecialFunctionGenerator_h
19 #define __itkSpecialFunctionGenerator_h
20 
21 #include "utlCore.h"
22 #include "utlNDArray.h"
23 
24 
25 namespace utl
26 {
34 inline double
35 Hyperg1F1(double a, double b, double x);
36 
37 
39 template < class T >
40 T
41 Lagurre ( const int n, const double a, const T x );
42 
43 
51 inline double
52 Gamma(const double x);
53 
54 inline double
55 GammaLower (const double s, const double x);
56 
66 inline double
67 BesselJa(const double a, const double x);
68 
69 inline double
70 BesselJInteger(const int n, const double x);
71 
72 inline double
73 BesselJIntegerPrime(const int n, const double x);
74 
76 // template<class T>
77 // NDArray<T,1>
78 // GetRotatedSHCoefficients(const NDArray<T,1>& shInput, const NDArray<T,2>& rotationMatrix);
79 
83 template < class T >
84 utl_shared_ptr<utl::NDArray<T,1> >
85 GetSymmetricTensorSHCoef(const T b, const T e1, const T e2, const int lMax, const T theta=0, const T phi=0 );
86 
92 template < class T >
93 std::vector< std::vector<T> >
94 GetSymmetricTensorSHCoefDerivative(const T b, const T e1, const T e2, const int lMax, const T theta=0, const T phi=0 );
95 
101 inline double
102 GetExpProductLegendreCoef(const double a, const double b, const int l );
103 
108 template < class T >
109 utl_shared_ptr<utl::NDArray<T,1> >
110 ComputeDWISHCoefficientsForGPDCylinder (const T radius, const T diffusivity, const T deltaBig, const T deltaSmall, const T qq, const int lMax , const T theta=0, const T phi=0);
111 
112 
119 template < class T >
120 double
122 
128 inline double
129 ComputeOrientationalOrderFromSymmetricTensor(const double e1, const double e2, const double phi=0);
130 
133 }
134 
135 
136 #if !defined(__itkSpecialFunctionGenerator_hxx)
138 #endif
139 
140 #endif
NDArray<T,1> is a vector class which uses blas mkl.
Definition: utlVector.h:36
double BesselJIntegerPrime(const int n, const double x)
double ComputeOrientationalOrderFromSymmetricTensor(const double e1, const double e2, const double phi=0)
double GammaLower(const double s, const double x)
double GetExpProductLegendreCoef(const double a, const double b, const int l)
double Hyperg1F1(double a, double b, double x)
Definition: utl.h:90
T Lagurre(const int n, const double a, const T x)
utl_shared_ptr< utl::NDArray< T, 1 > > ComputeDWISHCoefficientsForGPDCylinder(const T radius, const T diffusivity, const T deltaBig, const T deltaSmall, const T qq, const int lMax, const T theta=0, const T phi=0)
double BesselJInteger(const int n, const double x)
std::vector< std::vector< T > > GetSymmetricTensorSHCoefDerivative(const T b, const T e1, const T e2, const int lMax, const T theta=0, const T phi=0)
get the derivatives of the SH coefficients with respect to (e1, e2) in the symmetric tensor with eige...
double ComputeOrientationalOrderFromSHCoefficients(const utl::NDArray< T, 1 > &shCoef, const utl::NDArray< T, 1 > &axis)
double BesselJa(const double a, const double x)
bessel_Ja : Regular Cylindrical Bessel Function
utl_shared_ptr< utl::NDArray< T, 1 > > GetSymmetricTensorSHCoef(const T b, const T e1, const T e2, const int lMax, const T theta=0, const T phi=0)
get the SH coefficients from the symmetric tensor with eigenvalues (e1,e2,e2), e1>e2, and (theta,phi) is the angular direction of the e1 axis.
double Gamma(const double x)
gamma function.