65 template <
typename value_type,
class VectorType>
66 Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>
init_amp_transform (
const VectorType& els,
const size_t lmax)
68 Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic> ZSHT;
70 Eigen::Matrix<value_type, Eigen::Dynamic, 1, 0, 64> AL (lmax+1);
71 for (
size_t i = 0; i != size_t(els.size()); i++) {
73 for (
size_t l = 0; l <= lmax; l += 2)
74 ZSHT (i,
index(l)) = AL[l];
86 template <
typename value_type,
class VectorType>
87 Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic>
init_deriv_transform (
const VectorType& els,
const size_t lmax)
89 Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic> dZSHdelT;
90 dZSHdelT.resize (els.size(),
ZSH::NforL (lmax));
91 Eigen::Matrix<value_type, Eigen::Dynamic, 1, 0, 64> AL (lmax+1);
92 for (
size_t i = 0; i != size_t(els.size()); i++) {
94 dZSHdelT (i,
index(0)) = 0.0;
95 for (
size_t l = 2; l <= lmax; l += 2)
104 template <
typename ValueType>
107 using matrix_type = Eigen::Matrix<ValueType, Eigen::Dynamic, Eigen::Dynamic>;
109 template <
class MatrixType>
110 Transform (
const MatrixType& dirs,
const size_t lmax) :
114 template <
class VectorType1,
class VectorType2>
115 void A2ZSH (VectorType1& zsh,
const VectorType2& amplitudes)
const {
116 zsh.noalias() =
iZSHT * amplitudes;
118 template <
class VectorType1,
class VectorType2>
119 void ZSH2A (VectorType1& amplitudes,
const VectorType2& zsh)
const {
120 amplitudes.noalias() =
ZSHT * zsh;
123 size_t n_ZSH ()
const {
126 size_t n_amp ()
const {
143 template <
class VectorType>
144 inline typename VectorType::Scalar
value (
145 const VectorType& coefs,
146 typename VectorType::Scalar elevation,
149 using value_type =
typename VectorType::Scalar;
150 Eigen::Matrix<value_type, Eigen::Dynamic, 1, 0, 64> AL (lmax+1);
153 for (
size_t l = 0; l <= lmax; l += 2)
154 amplitude += AL[l] * coefs[
index(l)];
160 template <
class VectorType>
162 const VectorType& coefs,
163 const typename VectorType::Scalar elevation,
166 using value_type =
typename VectorType::Scalar;
167 Eigen::Matrix<value_type, Eigen::Dynamic, 1, 0, 64> AL (lmax+1);
170 for (
size_t l = 2; l <= lmax; l += 2)
177 template <
class VectorType1,
class VectorType2>
178 inline VectorType1&
ZSH2SH (VectorType1& sh,
const VectorType2& zsh)
180 const size_t lmax =
LforN (zsh.size());
182 for (
size_t i = 0; i != size_t(sh.size()); ++i)
184 for (
size_t l = 0; l <= lmax; l+=2)
189 template <
class VectorType>
190 inline Eigen::Matrix<typename VectorType::Scalar, Eigen::Dynamic, 1>
ZSH2SH (
const VectorType& zsh)
192 Eigen::Matrix<typename VectorType::Scalar, Eigen::Dynamic, 1> sh;
200 template <
class VectorType1,
class VectorType2>
201 inline VectorType1&
SH2ZSH (VectorType1& zsh,
const VectorType2& sh)
204 zsh.resize (
NforL (lmax));
205 for (
size_t l = 0; l <= lmax; l+=2)
210 template <
class VectorType>
211 inline Eigen::Matrix<typename VectorType::Scalar, Eigen::Dynamic, 1>
SH2ZSH (
const VectorType& sh)
213 Eigen::Matrix<typename VectorType::Scalar, Eigen::Dynamic, 1> zsh;
220 template <
class VectorType1,
class VectorType2>
221 inline VectorType1&
ZSH2RH (VectorType1& rh,
const VectorType2& zsh)
223 using value_type =
typename VectorType2::Scalar;
224 rh.resize (zsh.size());
225 const size_t lmax =
LforN (zsh.size());
226 Eigen::Matrix<value_type,Eigen::Dynamic,1,0,64> AL (lmax+1);
228 for (
size_t l = 0; l <= lmax; l += 2)
233 template <
class VectorType>
234 inline Eigen::Matrix<typename VectorType::Scalar,Eigen::Dynamic,1>
ZSH2RH (
const VectorType& zsh)
236 Eigen::Matrix<typename VectorType::Scalar,Eigen::Dynamic,1> rh (zsh.size());
246 template <
class VectorType1,
class VectorType2>
247 inline VectorType1&
zsconv (VectorType1& zsh,
const VectorType2& RH)
249 assert (zsh.size() >= RH.size());
250 for (
size_t i = 0; i != size_t(RH.size()); ++i)
259 template <
class VectorType1,
class VectorType2,
class VectorType3>
260 inline VectorType1&
zsconv (VectorType1& C,
const VectorType2& RH,
const VectorType3& zsh)
262 assert (zsh.size() >= RH.size());
263 C.resize (RH.size());
264 for (
size_t i = 0; i != size_t(RH.size()); ++i)
265 C[i] = zsh[i] * RH[i];
273 template <
class VectorType>
277 default_type ev1 = ADC* (1.0+2.0*a), ev2 = ADC* (1.0-a);
279 Eigen::VectorXd sigs (precision);
280 Eigen::MatrixXd ZSHT (precision, lmax/2+1);
281 Eigen::Matrix<default_type,Eigen::Dynamic,1,0,64> AL;
283 for (
size_t i = 0; i < precision; i++) {
285 sigs[i] = exp (-bvalue * (ev1*std::cos (
el)*std::cos (
el) + ev2*std::sin (
el)*std::sin (
el)));
287 for (
size_t l = 0; l <= lmax; l+=2)
288 ZSHT (i,
index(l)) = AL[l];
291 return (zsh =
pinv(ZSHT) * sigs);
Eigen::Matrix< typename MatrixType::Scalar, Eigen::Dynamic, Eigen::Dynamic > pinv(const MatrixType &M)
return Moore-Penrose pseudo-inverse of M
size_t index(int l, int m)
compute the index for coefficient (l,m)
size_t LforN(int N)
returns the largest lmax given N parameters
size_t NforL(int lmax)
the number of (even-degree) coefficients for the given value of lmax
VectorType1 & ZSH2RH(VectorType1 &rh, const VectorType2 &zsh)
VectorType1 & ZSH2SH(VectorType1 &sh, const VectorType2 &zsh)
VectorType::Scalar derivative(const VectorType &coefs, const typename VectorType::Scalar elevation, const size_t lmax)
size_t index(int l)
compute the index for coefficient l
VectorType & FA2ZSH(VectorType &zsh, default_type FA, default_type ADC, default_type bvalue, const size_t lmax, const size_t precision=100)
compute ZSH coefficients corresponding to specified tensor
size_t NforL(int lmax)
the number of (even-degree) coefficients for the given value of lmax
size_t LforN(int N)
returns the largest lmax given N parameters
VectorType1 & SH2ZSH(VectorType1 &zsh, const VectorType2 &sh)
VectorType1 & zsconv(VectorType1 &zsh, const VectorType2 &RH)
perform zonal spherical convolution, in place
Eigen::Matrix< value_type, Eigen::Dynamic, Eigen::Dynamic > init_deriv_transform(const VectorType &els, const size_t lmax)
form the ZSH->derivatives matrix for a set of elevation angles
VectorType::Scalar value(const VectorType &coefs, typename VectorType::Scalar elevation, const size_t lmax)
Eigen::Matrix< value_type, Eigen::Dynamic, Eigen::Dynamic > init_amp_transform(const VectorType &els, const size_t lmax)
form the ZSH->amplitudes matrix for a set of elevation angles
ValueType Plm_sph(const int l, const int m, const ValueType x)
MR::default_type value_type
Eigen::Matrix< value_type, Eigen::Dynamic, Eigen::Dynamic > matrix_type
double default_type
the default type used throughout MRtrix