Developer documentation
Version 3.0.3-105-gd3941f44
SH.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_SH_h__
18#define __math_SH_h__
19
20#include "math/legendre.h"
21#include "math/least_squares.h"
22
23#define MAX_DIR_CHANGE 0.2
24#define ANGLE_TOLERANCE 1e-4
25
26namespace MR
27{
28 namespace Math
29 {
30 namespace SH
31 {
32
41
43 extern const char* encoding_description;
44
46 inline size_t NforL (int lmax)
47 {
48 return (lmax+1) * (lmax+2) /2;
49 }
50
52 inline size_t index (int l, int m)
53 {
54 return l * (l+1) /2 + m;
55 }
56
58 inline size_t NforL_mpos (int lmax)
59 {
60 return (lmax/2+1) * (lmax/2+1);
61 }
62
64 inline size_t index_mpos (int l, int m)
65 {
66 return l*l/4 + m;
67 }
68
70 inline size_t LforN (int N)
71 {
72 return N ? 2 * std::floor<size_t> ( (std::sqrt (float (1+8*N))-3.0) /4.0) : 0;
73 }
74
75
76
78
81 template <class MatrixType>
82 Eigen::Matrix<typename MatrixType::Scalar,Eigen::Dynamic, Eigen::Dynamic> init_transform (const MatrixType& dirs, const int lmax)
83 {
84 using namespace Eigen;
85 using value_type = typename MatrixType::Scalar;
86 if (dirs.cols() != 2)
87 throw Exception ("direction matrix should have 2 columns: [ azimuth elevation ]");
88 Matrix<value_type,Dynamic,Dynamic> SHT (dirs.rows(), NforL (lmax));
89 Matrix<value_type,Dynamic,1,0,64> AL (lmax+1);
90 for (ssize_t i = 0; i < dirs.rows(); i++) {
91 const value_type z = std::cos (dirs (i,1));
92 Legendre::Plm_sph (AL, lmax, 0, z);
93 for (int l = 0; l <= lmax; l+=2)
94 SHT (i,index (l,0)) = AL[l];
95 for (int m = 1; m <= lmax; m++) {
96 Legendre::Plm_sph (AL, lmax, m, z);
97 for (int l = ( (m&1) ? m+1 : m); l <= lmax; l+=2) {
98 SHT(i, index(l, m)) = Math::sqrt2 * AL[l]*std::cos (m*dirs (i,0));
99 SHT(i, index(l,-m)) = Math::sqrt2 * AL[l]*std::sin (m*dirs (i,0));
100 }
101 }
102 }
103 return SHT;
104 }
105
107
110 template <class MatrixType>
111 Eigen::Matrix<typename MatrixType::Scalar,Eigen::Dynamic, Eigen::Dynamic> init_transform_cart (const MatrixType& dirs, const int lmax)
112 {
113 using namespace Eigen;
114 using value_type = typename MatrixType::Scalar;
115 if (dirs.cols() != 3)
116 throw Exception ("direction matrix should have 3 columns: [ x y z ]");
117 Matrix<value_type,Dynamic,Dynamic> SHT (dirs.rows(), NforL (lmax));
118 Matrix<value_type,Dynamic,1,0,64> AL (lmax+1);
119 for (ssize_t i = 0; i < dirs.rows(); i++) {
120 value_type z = dirs (i,2);
121 value_type rxy = std::hypot(dirs(i,0), dirs(i,1));
122 value_type cp = (rxy) ? dirs(i,0)/rxy : 1.0;
123 value_type sp = (rxy) ? dirs(i,1)/rxy : 0.0;
124 Legendre::Plm_sph (AL, lmax, 0, z);
125 for (int l = 0; l <= lmax; l+=2)
126 SHT (i,index (l,0)) = AL[l];
127 value_type c0 (1.0), s0 (0.0);
128 for (int m = 1; m <= lmax; m++) {
129 Legendre::Plm_sph (AL, lmax, m, z);
130 value_type c = c0 * cp - s0 * sp;
131 value_type s = s0 * cp + c0 * sp;
132 for (int l = ( (m&1) ? m+1 : m); l <= lmax; l+=2) {
133 SHT(i, index(l, m)) = Math::sqrt2 * AL[l] * c;
134 SHT(i, index(l,-m)) = Math::sqrt2 * AL[l] * s;
135 }
136 c0 = c;
137 s0 = s;
138 }
139 }
140 return SHT;
141 }
142
143
144
146 template <class MatrixType, class VectorType>
147 inline void scale_degrees_forward (MatrixType& SH2amp_mapping, const VectorType& coefs)
148 {
149 ssize_t l = 0, nl = 1;
150 for (ssize_t col = 0; col < SH2amp_mapping.cols(); ++col) {
151 if (col >= nl) {
152 l++;
153 nl = NforL (2*l);
154 }
155 SH2amp_mapping.col(col) *= coefs[l];
156 }
157 }
158
159
160
162 template <typename MatrixType, class VectorType>
163 inline void scale_degrees_inverse (MatrixType& amp2SH_mapping, const VectorType& coefs)
164 {
165 ssize_t l = 0, nl = 1;
166 for (ssize_t row = 0; row < amp2SH_mapping.rows(); ++row) {
167 if (row >= nl) {
168 l++;
169 nl = NforL (2*l);
170 }
171 amp2SH_mapping.row(row) *= coefs[l];
172 }
173 }
174
175
176
177
179 template <typename VectorType>
180 inline Eigen::Matrix<typename VectorType::Scalar,Eigen::Dynamic,1> invert (const VectorType& coefs)
181 {
182 Eigen::Matrix<typename VectorType::Scalar,Eigen::Dynamic,1> ret (coefs.size());
183 for (size_t n = 0; n < coefs.size(); ++n)
184 ret[n] = ( coefs[n] ? 1.0 / coefs[n] : 0.0 );
185 return ret;
186 }
187
188 template <typename ValueType>
190 public:
191 using matrix_type = Eigen::Matrix<ValueType,Eigen::Dynamic,Eigen::Dynamic>;
192
193 template <class MatrixType>
194 Transform (const MatrixType& dirs, int lmax) :
195 SHT (init_transform (dirs, lmax)),
196 iSHT (pinv (SHT)) { }
197
198 template <class VectorType>
199 void set_filter (const VectorType& filter) {
200 scale_degrees_forward (SHT, invert (filter));
201 scale_degrees_inverse (iSHT, filter);
202 }
203 template <class VectorType1, class VectorType2>
204 void A2SH (VectorType1& sh, const VectorType2& amplitudes) const {
205 sh.noalias() = iSHT * amplitudes;
206 }
207 template <class VectorType1, class VectorType2>
208 void SH2A (VectorType1& amplitudes, const VectorType2& sh) const {
209 amplitudes.noalias() = SHT * sh;
210 }
211
212 size_t n_SH () const {
213 return SHT.cols();
214 }
215 size_t n_amp () const {
216 return SHT.rows();
217 }
218
219 const matrix_type& mat_A2SH () const {
220 return iSHT;
221 }
222 const matrix_type& mat_SH2A () const {
223 return SHT;
224 }
225
226 protected:
228 };
229
230
231
232 template <class VectorType>
233 inline typename VectorType::Scalar value (const VectorType& coefs,
234 typename VectorType::Scalar cos_elevation,
235 typename VectorType::Scalar cos_azimuth,
236 typename VectorType::Scalar sin_azimuth,
237 int lmax)
238 {
239 using value_type = typename VectorType::Scalar;
240 value_type amplitude = 0.0;
241 Eigen::Matrix<value_type,Eigen::Dynamic,1,0,64> AL (lmax+1);
242 Legendre::Plm_sph (AL, lmax, 0, cos_elevation);
243 for (int l = 0; l <= lmax; l+=2)
244 amplitude += AL[l] * coefs[index (l,0)];
245 value_type c0 (1.0), s0 (0.0);
246 for (int m = 1; m <= lmax; m++) {
247 Legendre::Plm_sph (AL, lmax, m, cos_elevation);
248 value_type c = c0 * cos_azimuth - s0 * sin_azimuth; // std::cos(m*azimuth)
249 value_type s = s0 * cos_azimuth + c0 * sin_azimuth; // std::sin(m*azimuth)
250 for (int l = ( (m&1) ? m+1 : m); l <= lmax; l+=2)
251 amplitude += AL[l] * Math::sqrt2 * (c * coefs[index (l,m)] + s * coefs[index (l,-m)]);
252 c0 = c;
253 s0 = s;
254 }
255 return amplitude;
256 }
257
258 template <class VectorType>
259 inline typename VectorType::Scalar value (const VectorType& coefs,
260 typename VectorType::Scalar cos_elevation,
261 typename VectorType::Scalar azimuth,
262 int lmax)
263 {
264 return value (coefs, cos_elevation, std::cos(azimuth), std::sin(azimuth), lmax);
265 }
266
267 template <class VectorType1, class VectorType2>
268 inline typename VectorType1::Scalar value (const VectorType1& coefs, const VectorType2& unit_dir, int lmax)
269 {
270 using value_type = typename VectorType1::Scalar;
271 value_type rxy = std::sqrt ( pow2(unit_dir[1]) + pow2(unit_dir[0]) );
272 value_type cp = (rxy) ? unit_dir[0]/rxy : 1.0;
273 value_type sp = (rxy) ? unit_dir[1]/rxy : 0.0;
274 return value (coefs, unit_dir[2], cp, sp, lmax);
275 }
276
277
278 template <class VectorType1, class VectorType2>
279 inline VectorType1& delta (VectorType1& delta_vec, const VectorType2& unit_dir, int lmax)
280 {
281 using value_type = typename VectorType1::Scalar;
282 delta_vec.resize (NforL (lmax));
283 value_type rxy = std::sqrt ( pow2(unit_dir[1]) + pow2(unit_dir[0]) );
284 value_type cp = (rxy) ? unit_dir[0]/rxy : 1.0;
285 value_type sp = (rxy) ? unit_dir[1]/rxy : 0.0;
286 Eigen::Matrix<value_type,Eigen::Dynamic,1,0,64> AL (lmax+1);
287 Legendre::Plm_sph (AL, lmax, 0, unit_dir[2]);
288 for (int l = 0; l <= lmax; l+=2)
289 delta_vec[index (l,0)] = AL[l];
290 value_type c0 (1.0), s0 (0.0);
291 for (int m = 1; m <= lmax; m++) {
292 Legendre::Plm_sph (AL, lmax, m, unit_dir[2]);
293 value_type c = c0 * cp - s0 * sp;
294 value_type s = s0 * cp + c0 * sp;
295 for (int l = ( (m&1) ? m+1 : m); l <= lmax; l+=2) {
296 delta_vec[index (l,m)] = AL[l] * Math::sqrt2 * c;
297 delta_vec[index (l,-m)] = AL[l] * Math::sqrt2 * s;
298 }
299 c0 = c;
300 s0 = s;
301 }
302 return delta_vec;
303 }
304
305
306
307 template <class VectorType1, class VectorType2>
308 inline VectorType1& SH2RH (VectorType1& RH, const VectorType2& sh)
309 {
310 using value_type = typename VectorType2::Scalar;
311 RH.resize (sh.size());
312 int lmax = 2*sh.size() +1;
313 Eigen::Matrix<value_type,Eigen::Dynamic,1,0,64> AL (lmax+1);
314 Legendre::Plm_sph (AL, lmax, 0, 1.0);
315 for (ssize_t l = 0; l < sh.size(); l++)
316 RH[l] = sh[l]/ AL[2*l];
317 return RH;
318 }
319
320 template <class VectorType>
321 inline Eigen::Matrix<typename VectorType::Scalar,Eigen::Dynamic,1> SH2RH (const VectorType& sh)
322 {
323 Eigen::Matrix<typename VectorType::Scalar,Eigen::Dynamic,1> RH (sh.size());
324 SH2RH (RH, sh);
325 return RH;
326 }
327
328
330
332 template <class VectorType1, class VectorType2>
333 inline VectorType1& sconv (VectorType1& sh, const VectorType2& RH)
334 {
335 assert (sh.size() >= ssize_t (NforL (2* (RH.size()-1))));
336 for (ssize_t i = 0; i < RH.size(); ++i) {
337 int l = 2*i;
338 for (int m = -l; m <= l; ++m)
339 sh[index (l,m)] *= RH[i];
340 }
341 return sh;
342 }
343
344
346
348 template <class VectorType1, class VectorType2, class VectorType3>
349 inline VectorType1& sconv (VectorType1& C, const VectorType2& RH, const VectorType3& sh)
350 {
351 assert (sh.size() >= ssize_t (NforL (2* (RH.size()-1))));
352 C.resize (NforL (2* (RH.size()-1)));
353 for (ssize_t i = 0; i < RH.size(); ++i) {
354 int l = 2*i;
355 for (int m = -l; m <= l; ++m)
356 C[index (l,m)] = RH[i] * sh[index (l,m)];
357 }
358 return C;
359 }
360
361
363
366 template <class MatrixType1, class VectorType2>
367 inline MatrixType1& sconv_mat (MatrixType1& sh, const VectorType2& RH)
368 {
369 assert (sh.cols() >= ssize_t (NforL (2* (RH.size()-1))));
370 for (ssize_t i = 0; i < RH.size(); ++i) {
371 int l = 2*i;
372 for (int m = -l; m <= l; ++m)
373 sh.col(index (l,m)) *= RH[i];
374 }
375 return sh;
376 }
377
378
379 namespace {
380 template <typename> struct __dummy { NOMEMALIGN using type = int; };
381 }
382
383
384
385
386
387
389 template <typename ValueType> class PrecomputedFraction
390 { NOMEMALIGN
391 public:
392 PrecomputedFraction () : f1 (0.0), f2 (0.0) { }
393 ValueType f1, f2;
395 };
396
397
399 template <typename ValueType> class PrecomputedAL
400 { NOMEMALIGN
401 public:
402 using value_type = ValueType;
403
404 PrecomputedAL () : lmax (0), ndir (0), nAL (0), inc (0.0) { }
405 PrecomputedAL (int up_to_lmax, int num_dir = 512) {
406 init (up_to_lmax, num_dir);
407 }
408
409 bool operator! () const {
410 return AL.empty();
411 }
412 operator bool () const {
413 return AL.size();
414 }
415
416 void init (int up_to_lmax, int num_dir = 512) {
417 lmax = up_to_lmax;
418 ndir = num_dir;
419 nAL = NforL_mpos (lmax);
420 inc = Math::pi / (ndir-1);
421 AL.resize (ndir*nAL);
422 Eigen::Matrix<value_type,Eigen::Dynamic,1,0,64> buf (lmax+1);
423
424 for (int n = 0; n < ndir; n++) {
425 typename vector<value_type>::iterator p = AL.begin() + n*nAL;
426 value_type cos_el = std::cos (n*inc);
427 for (int m = 0; m <= lmax; m++) {
428 Legendre::Plm_sph (buf, lmax, m, cos_el);
429 for (int l = ( (m&1) ?m+1:m); l <= lmax; l+=2)
430 p[index_mpos (l,m)] = buf[l];
431 }
432 }
433 }
434
435 void set (PrecomputedFraction<ValueType>& f, const ValueType elevation) const {
436 f.f2 = elevation / inc;
437 int i = int (f.f2);
438 if (i < 0) {
439 i = 0;
440 f.f1 = 1.0;
441 f.f2 = 0.0;
442 }
443 else if (i >= ndir-1) {
444 i = ndir-1;
445 f.f1 = 1.0;
446 f.f2 = 0.0;
447 }
448 else {
449 f.f2 -= i;
450 f.f1 = 1.0 - f.f2;
451 }
452
453 f.p1 = AL.begin() + i*nAL;
454 f.p2 = f.p1 + nAL;
455 }
456
457 ValueType get (const PrecomputedFraction<ValueType>& f, int i) const {
458 ValueType v = f.f1*f.p1[i];
459 if (f.f2) v += f.f2*f.p2[i];
460 return v;
461 }
462 ValueType get (const PrecomputedFraction<ValueType>& f, int l, int m) const {
463 return get (f, index_mpos (l,m));
464 }
465
466 void get (ValueType* dest, const PrecomputedFraction<ValueType>& f) const {
467 for (int l = 0; l <= lmax; l+=2) {
468 for (int m = 0; m <= l; m++) {
469 int i = index_mpos (l,m);
470 dest[i] = get (f,i);
471 }
472 }
473 }
474
475 template <class VectorType, class UnitVectorType>
476 ValueType value (const VectorType& val, const UnitVectorType& unit_dir) const {
478 set (f, std::acos (unit_dir[2]));
479 ValueType rxy = std::sqrt ( pow2(unit_dir[1]) + pow2(unit_dir[0]) );
480 ValueType cp = (rxy) ? unit_dir[0]/rxy : 1.0;
481 ValueType sp = (rxy) ? unit_dir[1]/rxy : 0.0;
482 ValueType v = 0.0;
483 for (int l = 0; l <= lmax; l+=2)
484 v += get (f,l,0) * val[index (l,0)];
485 ValueType c0 (1.0), s0 (0.0);
486 for (int m = 1; m <= lmax; m++) {
487 ValueType c = c0 * cp - s0 * sp;
488 ValueType s = s0 * cp + c0 * sp;
489 for (int l = ( (m&1) ? m+1 : m); l <= lmax; l+=2)
490 v += get (f,l,m) * Math::sqrt2 * (c * val[index (l,m)] + s * val[index (l,-m)]);
491 c0 = c;
492 s0 = s;
493 }
494 return v;
495 }
496
497 protected:
498 int lmax, ndir, nAL;
499 ValueType inc;
501 };
502
503
504
505
506
507
509
514 template <class VectorType, class UnitVectorType, class ValueType = float>
515 inline typename VectorType::Scalar get_peak (
516 const VectorType& sh,
517 int lmax,
518 UnitVectorType& unit_init_dir,
519 PrecomputedAL<typename VectorType::Scalar>* precomputer = nullptr)
520 {
521 using value_type = typename VectorType::Scalar;
522 assert (std::isfinite (unit_init_dir[0]));
523 for (int i = 0; i < 50; i++) {
524 value_type az = std::atan2 (unit_init_dir[1], unit_init_dir[0]);
525 value_type el = std::acos (unit_init_dir[2]);
526 value_type amplitude, dSH_del, dSH_daz, d2SH_del2, d2SH_deldaz, d2SH_daz2;
527 derivatives (sh, lmax, el, az, amplitude, dSH_del, dSH_daz, d2SH_del2, d2SH_deldaz, d2SH_daz2, precomputer);
528
529 value_type del = sqrt (dSH_del*dSH_del + dSH_daz*dSH_daz);
530 value_type daz = 0.0;
531 if (del != 0.0) {
532 daz = dSH_daz/del;
533 del = dSH_del/del;
534 }
535
536
537 value_type dSH_dt = daz*dSH_daz + del*dSH_del;
538 value_type d2SH_dt2 = daz*daz*d2SH_daz2 + 2.0*daz*del*d2SH_deldaz + del*del*d2SH_del2;
539 value_type dt = d2SH_dt2 ? (-dSH_dt / d2SH_dt2) : 0.0;
540
541 if (dt < 0.0) dt = -dt;
542 if (dt > MAX_DIR_CHANGE) dt = MAX_DIR_CHANGE;
543
544 del *= dt;
545 daz *= dt;
546
547 unit_init_dir[0] += del*std::cos (az) *std::cos (el) - daz*std::sin (az);
548 unit_init_dir[1] += del*std::sin (az) *std::cos (el) + daz*std::cos (az);
549 unit_init_dir[2] -= del*std::sin (el);
550 unit_init_dir.normalize();
551
552 if (dt < ANGLE_TOLERANCE)
553 return amplitude;
554 }
555
556 unit_init_dir = { NaN, NaN, NaN };
557 DEBUG ("failed to find SH peak!");
558 return NaN;
559 }
560
561
562
563
564
565
567
568 template <class VectorType>
569 inline void derivatives (
570 const VectorType& sh,
571 const int lmax,
572 const typename VectorType::Scalar elevation,
573 const typename VectorType::Scalar azimuth,
574 typename VectorType::Scalar& amplitude,
575 typename VectorType::Scalar& dSH_del,
576 typename VectorType::Scalar& dSH_daz,
577 typename VectorType::Scalar& d2SH_del2,
578 typename VectorType::Scalar& d2SH_deldaz,
579 typename VectorType::Scalar& d2SH_daz2,
581 {
582 using value_type = typename VectorType::Scalar;
583 value_type sel = std::sin (elevation);
584 value_type cel = std::cos (elevation);
585 bool atpole = sel < 1e-4;
586
587 dSH_del = dSH_daz = d2SH_del2 = d2SH_deldaz = d2SH_daz2 = 0.0;
588 VLA_MAX (AL, value_type, NforL_mpos (lmax), 64);
589
590 if (precomputer) {
592 precomputer->set (f, elevation);
593 precomputer->get (AL, f);
594 }
595 else {
596 Eigen::Matrix<value_type,Eigen::Dynamic,1,0,64> buf (lmax+1);
597 for (int m = 0; m <= lmax; m++) {
598 Legendre::Plm_sph (buf, lmax, m, cel);
599 for (int l = ( (m&1) ?m+1:m); l <= lmax; l+=2)
600 AL[index_mpos (l,m)] = buf[l];
601 }
602 }
603
604 amplitude = sh[0] * AL[0];
605 for (int l = 2; l <= (int) lmax; l+=2) {
606 const value_type& v (sh[index (l,0)]);
607 amplitude += v * AL[index_mpos (l,0)];
608 dSH_del += v * sqrt (value_type (l* (l+1))) * AL[index_mpos (l,1)];
609 d2SH_del2 += v * (sqrt (value_type (l* (l+1) * (l-1) * (l+2))) * AL[index_mpos (l,2)] - l* (l+1) * AL[index_mpos (l,0)]) /2.0;
610 }
611
612 for (int m = 1; m <= lmax; m++) {
613 value_type caz = Math::sqrt2 * std::cos (m*azimuth);
614 value_type saz = Math::sqrt2 * std::sin (m*azimuth);
615 for (int l = ( (m&1) ? m+1 : m); l <= lmax; l+=2) {
616 const value_type& vp (sh[index (l,m)]);
617 const value_type& vm (sh[index (l,-m)]);
618 amplitude += (vp*caz + vm*saz) * AL[index_mpos (l,m)];
619
620 value_type tmp = sqrt (value_type ( (l+m) * (l-m+1))) * AL[index_mpos (l,m-1)];
621 if (l > m) tmp -= sqrt (value_type ( (l-m) * (l+m+1))) * AL[index_mpos (l,m+1)];
622 tmp /= -2.0;
623 dSH_del += (vp*caz + vm*saz) * tmp;
624
625 value_type tmp2 = - ( (l+m) * (l-m+1) + (l-m) * (l+m+1)) * AL[index_mpos (l,m)];
626 if (m == 1) tmp2 -= sqrt (value_type ( (l+m) * (l-m+1) * (l+m-1) * (l-m+2))) * AL[index_mpos (l,1)];
627 else tmp2 += sqrt (value_type ( (l+m) * (l-m+1) * (l+m-1) * (l-m+2))) * AL[index_mpos (l,m-2)];
628 if (l > m+1) tmp2 += sqrt (value_type ( (l-m) * (l+m+1) * (l-m-1) * (l+m+2))) * AL[index_mpos (l,m+2)];
629 tmp2 /= 4.0;
630 d2SH_del2 += (vp*caz + vm*saz) * tmp2;
631
632 if (atpole) dSH_daz += (vm*caz - vp*saz) * tmp;
633 else {
634 d2SH_deldaz += m * (vm*caz - vp*saz) * tmp;
635 dSH_daz += m * (vm*caz - vp*saz) * AL[index_mpos (l,m)];
636 d2SH_daz2 -= (vp*caz + vm*saz) * m*m * AL[index_mpos (l,m)];
637 }
638
639 }
640 }
641
642 if (!atpole) {
643 dSH_daz /= sel;
644 d2SH_deldaz /= sel;
645 d2SH_daz2 /= sel*sel;
646 }
647 }
648
649
650
652 template <typename ValueType> class aPSF
653 { MEMALIGN(aPSF<ValueType>)
654 public:
655 aPSF (const size_t lmax) :
656 lmax (lmax),
657 RH (lmax/2 + 1)
658 {
659 switch (lmax) {
660 case 2:
661 RH[0] = ValueType(1.00000000);
662 RH[1] = ValueType(0.41939279);
663 break;
664 case 4:
665 RH[0] = ValueType(1.00000000);
666 RH[1] = ValueType(0.63608543);
667 RH[2] = ValueType(0.18487087);
668 break;
669 case 6:
670 RH[0] = ValueType(1.00000000);
671 RH[1] = ValueType(0.75490341);
672 RH[2] = ValueType(0.37126442);
673 RH[3] = ValueType(0.09614699);
674 break;
675 case 8:
676 RH[0] = ValueType(1.00000000);
677 RH[1] = ValueType(0.82384816);
678 RH[2] = ValueType(0.51261696);
679 RH[3] = ValueType(0.22440563);
680 RH[4] = ValueType(0.05593079);
681 break;
682 case 10:
683 RH[0] = ValueType(1.00000000);
684 RH[1] = ValueType(0.86725945);
685 RH[2] = ValueType(0.61519436);
686 RH[3] = ValueType(0.34570667);
687 RH[4] = ValueType(0.14300355);
688 RH[5] = ValueType(0.03548062);
689 break;
690 case 12:
691 RH[0] = ValueType(1.00000000);
692 RH[1] = ValueType(0.89737759);
693 RH[2] = ValueType(0.69278503);
694 RH[3] = ValueType(0.45249879);
695 RH[4] = ValueType(0.24169922);
696 RH[5] = ValueType(0.09826171);
697 RH[6] = ValueType(0.02502481);
698 break;
699 case 14:
700 RH[0] = ValueType(1.00000000);
701 RH[1] = ValueType(0.91717853);
702 RH[2] = ValueType(0.74685644);
703 RH[3] = ValueType(0.53467773);
704 RH[4] = ValueType(0.33031863);
705 RH[5] = ValueType(0.17013825);
706 RH[6] = ValueType(0.06810155);
707 RH[7] = ValueType(0.01754930);
708 break;
709 case 16:
710 RH[0] = ValueType(1.00000000);
711 RH[1] = ValueType(0.93261196);
712 RH[2] = ValueType(0.79064858);
713 RH[3] = ValueType(0.60562880);
714 RH[4] = ValueType(0.41454703);
715 RH[5] = ValueType(0.24880754);
716 RH[6] = ValueType(0.12661242);
717 RH[7] = ValueType(0.05106681);
718 RH[8] = ValueType(0.01365433);
719 break;
720 default:
721 throw Exception ("No aPSF RH data for lmax " + str(lmax));
722 }
723 }
724
725
726 template <class VectorType, class UnitVectorType>
727 VectorType& operator() (VectorType& sh, const UnitVectorType& dir) const
728 {
729 sh.resize (RH.size());
730 delta (sh, dir, lmax);
731 return sconv (sh, RH);
732 }
733
734 inline const Eigen::Matrix<ValueType,Eigen::Dynamic,1>& RH_coefs() const { return RH; }
735
736 private:
737 const size_t lmax;
738 Eigen::Matrix<ValueType,Eigen::Dynamic,1> RH;
739
740 };
741
742
744 template <class ImageType>
745 void check (const ImageType& H) {
746 if (H.ndim() < 4)
747 throw Exception ("image \"" + H.name() + "\" does not contain SH coefficients - not 4D");
748 size_t l = LforN (H.size(3));
749 if (l%2 || NforL (l) != size_t (H.size(3)))
750 throw Exception ("image \"" + H.name() + "\" does not contain SH coefficients - unexpected number of coefficients");
751 }
754 }
755 }
756}
757
758#endif
#define MAX_DIR_CHANGE
Definition: SH.h:23
#define ANGLE_TOLERANCE
Definition: SH.h:24
float el
Definition: calibrator.h:61
Precomputed Associated Legrendre Polynomials - used to speed up SH calculation.
Definition: SH.h:400
ValueType get(const PrecomputedFraction< ValueType > &f, int l, int m) const
Definition: SH.h:462
void get(ValueType *dest, const PrecomputedFraction< ValueType > &f) const
Definition: SH.h:466
void set(PrecomputedFraction< ValueType > &f, const ValueType elevation) const
Definition: SH.h:435
ValueType value_type
Definition: SH.h:402
bool operator!() const
Definition: SH.h:409
ValueType get(const PrecomputedFraction< ValueType > &f, int i) const
Definition: SH.h:457
PrecomputedAL(int up_to_lmax, int num_dir=512)
Definition: SH.h:405
void init(int up_to_lmax, int num_dir=512)
Definition: SH.h:416
ValueType value(const VectorType &val, const UnitVectorType &unit_dir) const
Definition: SH.h:476
vector< ValueType > AL
Definition: SH.h:500
used to speed up SH calculation
Definition: SH.h:390
vector< ValueType >::const_iterator p1
Definition: SH.h:394
vector< ValueType >::const_iterator p2
Definition: SH.h:394
matrix_type iSHT
Definition: SH.h:227
matrix_type SHT
Definition: SH.h:227
a class to hold the coefficients for an apodised point-spread function.
Definition: SH.h:653
#define DEBUG(msg)
Definition: exception.h:75
#define VLA_MAX(name, type, num, max)
Definition: types.h:119
constexpr T pow2(const T &v)
Definition: math.h:53
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 sqrt2
Definition: math.h:43
constexpr double e
Definition: math.h:39
constexpr double pi
Definition: math.h:40
void check(const ImageType &H)
convenience function to check if an input image can contain SH coefficients
Definition: SH.h:745
size_t index(int l, int m)
compute the index for coefficient (l,m)
Definition: SH.h:52
void scale_degrees_inverse(MatrixType &amp2SH_mapping, const VectorType &coefs)
scale the coefficients of each SH degree by the corresponding value in coefs
Definition: SH.h:163
size_t LforN(int N)
returns the largest lmax given N parameters
Definition: SH.h:70
void scale_degrees_forward(MatrixType &SH2amp_mapping, const VectorType &coefs)
scale the coefficients of each SH degree by the corresponding value in coefs
Definition: SH.h:147
Eigen::Matrix< typename MatrixType::Scalar, Eigen::Dynamic, Eigen::Dynamic > init_transform_cart(const MatrixType &dirs, const int lmax)
form the SH->amplitudes matrix
Definition: SH.h:111
void derivatives(const VectorType &sh, const int lmax, const typename VectorType::Scalar elevation, const typename VectorType::Scalar azimuth, typename VectorType::Scalar &amplitude, typename VectorType::Scalar &dSH_del, typename VectorType::Scalar &dSH_daz, typename VectorType::Scalar &d2SH_del2, typename VectorType::Scalar &d2SH_deldaz, typename VectorType::Scalar &d2SH_daz2, PrecomputedAL< typename VectorType::Scalar > *precomputer)
computes first and second order derivatives of SH series
Definition: SH.h:569
VectorType1 & sconv(VectorType1 &sh, const VectorType2 &RH)
perform spherical convolution, in place
Definition: SH.h:333
size_t index_mpos(int l, int m)
same as index(), but consider only non-negative orders m
Definition: SH.h:64
VectorType1 & delta(VectorType1 &delta_vec, const VectorType2 &unit_dir, int lmax)
Definition: SH.h:279
size_t NforL_mpos(int lmax)
same as NforL(), but consider only non-negative orders m
Definition: SH.h:58
size_t NforL(int lmax)
the number of (even-degree) coefficients for the given value of lmax
Definition: SH.h:46
MatrixType1 & sconv_mat(MatrixType1 &sh, const VectorType2 &RH)
perform spherical convolution, in place
Definition: SH.h:367
const char * encoding_description
a string containing a description of the SH storage convention
VectorType::Scalar get_peak(const VectorType &sh, int lmax, UnitVectorType &unit_init_dir, PrecomputedAL< typename VectorType::Scalar > *precomputer=nullptr)
estimate direction & amplitude of SH peak
Definition: SH.h:515
VectorType::Scalar value(const VectorType &coefs, typename VectorType::Scalar cos_elevation, typename VectorType::Scalar cos_azimuth, typename VectorType::Scalar sin_azimuth, int lmax)
Definition: SH.h:233
Eigen::Matrix< typename MatrixType::Scalar, Eigen::Dynamic, Eigen::Dynamic > init_transform(const MatrixType &dirs, const int lmax)
form the SH->amplitudes matrix
Definition: SH.h:82
VectorType1 & SH2RH(VectorType1 &RH, const VectorType2 &sh)
Definition: SH.h:308
Eigen::Matrix< typename VectorType::Scalar, Eigen::Dynamic, 1 > invert(const VectorType &coefs)
invert any non-zero coefficients in coefs
Definition: SH.h:180
#define NOMEMALIGN
Definition: memory.h:22
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
std::string str(const T &value, int precision=0)
Definition: mrtrix.h:247
constexpr default_type NaN
Definition: types.h:230
Eigen::MatrixXd H
#define MEMALIGN(...)
Definition: types.h:185