Developer documentation
Version 3.0.3-105-gd3941f44
ZSH.h
Go to the documentation of this file.
1/* Copyright (c) 2008-2022 the MRtrix3 contributors.
2 *
3 * This Source Code Form is subject to the terms of the Mozilla Public
4 * License, v. 2.0. If a copy of the MPL was not distributed with this
5 * file, You can obtain one at http://mozilla.org/MPL/2.0/.
6 *
7 * Covered Software is provided under this License on an "as is"
8 * basis, without warranty of any kind, either expressed, implied, or
9 * statutory, including, without limitation, warranties that the
10 * Covered Software is free of defects, merchantable, fit for a
11 * particular purpose or non-infringing.
12 * See the Mozilla Public License v. 2.0 for more details.
13 *
14 * For more details, see http://www.mrtrix.org/.
15 */
16
17#ifndef __math_ZSH_h__
18#define __math_ZSH_h__
19
20#include <Eigen/Dense>
21
22#include "math/legendre.h"
23#include "math/least_squares.h"
24#include "math/SH.h"
25
26namespace MR
27{
28 namespace Math
29 {
30 namespace ZSH
31 {
32
42 inline size_t NforL (int lmax)
43 {
44 return (1 + lmax/2);
45 }
46
48 inline size_t index (int l)
49 {
50 return (l/2);
51 }
52
54 inline size_t LforN (int N)
55 {
56 return (2 * (N-1));
57 }
58
59
60
62
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)
67 {
68 Eigen::Matrix<value_type, Eigen::Dynamic, Eigen::Dynamic> ZSHT;
69 ZSHT.resize (els.size(), ZSH::NforL (lmax));
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++) {
72 Legendre::Plm_sph (AL, lmax, 0, std::cos (els[i]));
73 for (size_t l = 0; l <= lmax; l += 2)
74 ZSHT (i,index(l)) = AL[l];
75 }
76 return ZSHT;
77 }
78
79
80
82
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)
88 {
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++) {
93 Legendre::Plm_sph (AL, lmax, 1, std::cos (els[i]));
94 dZSHdelT (i,index(0)) = 0.0;
95 for (size_t l = 2; l <= lmax; l += 2)
96 dZSHdelT (i,index(l)) = AL[l] * sqrt (value_type (l*(l+1)));
97 }
98 return dZSHdelT;
99 }
100
101
102
103
104 template <typename ValueType>
106 public:
107 using matrix_type = Eigen::Matrix<ValueType, Eigen::Dynamic, Eigen::Dynamic>;
108
109 template <class MatrixType>
110 Transform (const MatrixType& dirs, const size_t lmax) :
111 ZSHT (init_amp_transform (dirs.col(1), lmax)), // Elevation angles are second column of aximuth/elevation matrix
112 iZSHT (pinv (ZSHT)) { }
113
114 template <class VectorType1, class VectorType2>
115 void A2ZSH (VectorType1& zsh, const VectorType2& amplitudes) const {
116 zsh.noalias() = iZSHT * amplitudes;
117 }
118 template <class VectorType1, class VectorType2>
119 void ZSH2A (VectorType1& amplitudes, const VectorType2& zsh) const {
120 amplitudes.noalias() = ZSHT * zsh;
121 }
122
123 size_t n_ZSH () const {
124 return ZSHT.cols();
125 }
126 size_t n_amp () const {
127 return ZSHT.rows();
128 }
129
130 const matrix_type& mat_A2ZSH () const {
131 return iZSHT;
132 }
133 const matrix_type& mat_ZSH2A () const {
134 return ZSHT;
135 }
136
137 protected:
139 };
140
141
142
143 template <class VectorType>
144 inline typename VectorType::Scalar value (
145 const VectorType& coefs,
146 typename VectorType::Scalar elevation,
147 const size_t lmax)
148 {
149 using value_type = typename VectorType::Scalar;
150 Eigen::Matrix<value_type, Eigen::Dynamic, 1, 0, 64> AL (lmax+1);
151 Legendre::Plm_sph (AL, lmax, 0, std::cos(elevation));
152 value_type amplitude = 0.0;
153 for (size_t l = 0; l <= lmax; l += 2)
154 amplitude += AL[l] * coefs[index(l)];
155 return amplitude;
156 }
157
158
159
160 template <class VectorType>
161 inline typename VectorType::Scalar derivative (
162 const VectorType& coefs,
163 const typename VectorType::Scalar elevation,
164 const size_t lmax)
165 {
166 using value_type = typename VectorType::Scalar;
167 Eigen::Matrix<value_type, Eigen::Dynamic, 1, 0, 64> AL (lmax+1);
168 Legendre::Plm_sph (AL, lmax, 1, std::cos (elevation));
169 value_type dZSH_del = 0.0;
170 for (size_t l = 2; l <= lmax; l += 2)
171 dZSH_del += AL[l] * coefs[index(l)] * sqrt (value_type (l*(l+1)));
172 return dZSH_del;
173 }
174
175
176
177 template <class VectorType1, class VectorType2>
178 inline VectorType1& ZSH2SH (VectorType1& sh, const VectorType2& zsh)
179 {
180 const size_t lmax = LforN (zsh.size());
181 sh.resize (Math::SH::NforL (lmax));
182 for (size_t i = 0; i != size_t(sh.size()); ++i)
183 sh[i] = 0.0;
184 for (size_t l = 0; l <= lmax; l+=2)
185 sh[Math::SH::index(l,0)] = zsh[index(l)];
186 return sh;
187 }
188
189 template <class VectorType>
190 inline Eigen::Matrix<typename VectorType::Scalar, Eigen::Dynamic, 1> ZSH2SH (const VectorType& zsh)
191 {
192 Eigen::Matrix<typename VectorType::Scalar, Eigen::Dynamic, 1> sh;
193 ZSH2SH (sh, zsh);
194 return sh;
195 }
196
197
198
199
200 template <class VectorType1, class VectorType2>
201 inline VectorType1& SH2ZSH (VectorType1& zsh, const VectorType2& sh)
202 {
203 const size_t lmax = Math::SH::LforN (sh.size());
204 zsh.resize (NforL (lmax));
205 for (size_t l = 0; l <= lmax; l+=2)
206 zsh[index(l)] = sh[Math::SH::index(l,0)];
207 return zsh;
208 }
209
210 template <class VectorType>
211 inline Eigen::Matrix<typename VectorType::Scalar, Eigen::Dynamic, 1> SH2ZSH (const VectorType& sh)
212 {
213 Eigen::Matrix<typename VectorType::Scalar, Eigen::Dynamic, 1> zsh;
214 SH2ZSH (zsh, sh);
215 return zsh;
216 }
217
218
219
220 template <class VectorType1, class VectorType2>
221 inline VectorType1& ZSH2RH (VectorType1& rh, const VectorType2& zsh)
222 {
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);
227 Legendre::Plm_sph (AL, lmax, 0, 1.0);
228 for (size_t l = 0; l <= lmax; l += 2)
229 rh[index(l)] = zsh[index(l)] / AL[l];
230 return rh;
231 }
232
233 template <class VectorType>
234 inline Eigen::Matrix<typename VectorType::Scalar,Eigen::Dynamic,1> ZSH2RH (const VectorType& zsh)
235 {
236 Eigen::Matrix<typename VectorType::Scalar,Eigen::Dynamic,1> rh (zsh.size());
237 ZSH2RH (rh, zsh);
238 return rh;
239 }
240
241
242
244
246 template <class VectorType1, class VectorType2>
247 inline VectorType1& zsconv (VectorType1& zsh, const VectorType2& RH)
248 {
249 assert (zsh.size() >= RH.size());
250 for (size_t i = 0; i != size_t(RH.size()); ++i)
251 zsh[i] *= RH[i];
252 return zsh;
253 }
254
255
257
259 template <class VectorType1, class VectorType2, class VectorType3>
260 inline VectorType1& zsconv (VectorType1& C, const VectorType2& RH, const VectorType3& zsh)
261 {
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];
266 return C;
267 }
268
269
270
271
273 template <class VectorType>
274 inline VectorType& FA2ZSH (VectorType& zsh, default_type FA, default_type ADC, default_type bvalue, const size_t lmax, const size_t precision = 100)
275 {
276 default_type a = FA/sqrt (3.0 - 2.0*FA*FA);
277 default_type ev1 = ADC* (1.0+2.0*a), ev2 = ADC* (1.0-a);
278
279 Eigen::VectorXd sigs (precision);
280 Eigen::MatrixXd ZSHT (precision, lmax/2+1);
281 Eigen::Matrix<default_type,Eigen::Dynamic,1,0,64> AL;
282
283 for (size_t i = 0; i < precision; i++) {
284 default_type el = i*Math::pi / (2.0*(precision-1));
285 sigs[i] = exp (-bvalue * (ev1*std::cos (el)*std::cos (el) + ev2*std::sin (el)*std::sin (el)));
286 Legendre::Plm_sph (AL, lmax, 0, std::cos (el));
287 for (size_t l = 0; l <= lmax; l+=2)
288 ZSHT (i,index(l)) = AL[l];
289 }
290
291 return (zsh = pinv(ZSHT) * sigs);
292 }
293
294
295
296 }
297 }
298}
299
300#endif
float el
Definition: calibrator.h:61
Eigen::Matrix< typename MatrixType::Scalar, Eigen::Dynamic, Eigen::Dynamic > pinv(const MatrixType &M)
return Moore-Penrose pseudo-inverse of M
Definition: least_squares.h:39
constexpr double pi
Definition: math.h:40
size_t index(int l, int m)
compute the index for coefficient (l,m)
Definition: SH.h:52
size_t LforN(int N)
returns the largest lmax given N parameters
Definition: SH.h:70
size_t NforL(int lmax)
the number of (even-degree) coefficients for the given value of lmax
Definition: SH.h:46
VectorType1 & ZSH2RH(VectorType1 &rh, const VectorType2 &zsh)
Definition: ZSH.h:221
VectorType1 & ZSH2SH(VectorType1 &sh, const VectorType2 &zsh)
Definition: ZSH.h:178
VectorType::Scalar derivative(const VectorType &coefs, const typename VectorType::Scalar elevation, const size_t lmax)
Definition: ZSH.h:161
matrix_type ZSHT
Definition: ZSH.h:138
size_t index(int l)
compute the index for coefficient l
Definition: ZSH.h:48
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
Definition: ZSH.h:274
matrix_type iZSHT
Definition: ZSH.h:138
size_t NforL(int lmax)
the number of (even-degree) coefficients for the given value of lmax
Definition: ZSH.h:42
size_t LforN(int N)
returns the largest lmax given N parameters
Definition: ZSH.h:54
VectorType1 & SH2ZSH(VectorType1 &zsh, const VectorType2 &sh)
Definition: ZSH.h:201
VectorType1 & zsconv(VectorType1 &zsh, const VectorType2 &RH)
perform zonal spherical convolution, in place
Definition: ZSH.h:247
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
Definition: ZSH.h:87
VectorType::Scalar value(const VectorType &coefs, typename VectorType::Scalar elevation, const size_t lmax)
Definition: ZSH.h:144
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
Definition: ZSH.h:66
ValueType Plm_sph(const int l, const int m, const ValueType x)
Definition: legendre.h:73
MR::default_type value_type
Definition: typedefs.h:33
Eigen::Matrix< value_type, Eigen::Dynamic, Eigen::Dynamic > matrix_type
Definition: typedefs.h:34
Definition: base.h:24
double default_type
the default type used throughout MRtrix
Definition: types.h:228
#define MEMALIGN(...)
Definition: types.h:185