File: itkSymmetricEigenAnalysis.h

package info (click to toggle)
insighttoolkit5 5.4.3-5
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 704,384 kB
  • sloc: cpp: 783,592; ansic: 628,724; xml: 44,704; fortran: 34,250; python: 22,874; sh: 4,078; pascal: 2,636; lisp: 2,158; makefile: 464; yacc: 328; asm: 205; perl: 203; lex: 146; tcl: 132; javascript: 98; csh: 81
file content (1104 lines) | stat: -rw-r--r-- 41,640 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
/*=========================================================================
 *
 *  Copyright NumFOCUS
 *
 *  Licensed under the Apache License, Version 2.0 (the "License");
 *  you may not use this file except in compliance with the License.
 *  You may obtain a copy of the License at
 *
 *         https://www.apache.org/licenses/LICENSE-2.0.txt
 *
 *  Unless required by applicable law or agreed to in writing, software
 *  distributed under the License is distributed on an "AS IS" BASIS,
 *  WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
 *  See the License for the specific language governing permissions and
 *  limitations under the License.
 *
 *=========================================================================*/
#ifndef itkSymmetricEigenAnalysis_h
#define itkSymmetricEigenAnalysis_h

#include "itkMacro.h"
#include "itk_eigen.h"
#include ITK_EIGEN(Eigenvalues)
#include <numeric>
#include <vector>
// For GetPointerToMatrixData
#include "vnl/vnl_matrix.h"
#include "vnl/vnl_matrix_fixed.h"
#include "itkMatrix.h"

namespace itk
{
namespace detail
{
/* Helper functions returning pointer to matrix data for different types.  */
template <typename TValueType, unsigned int VRows, unsigned int VColumns>
const TValueType *
GetPointerToMatrixData(const vnl_matrix_fixed<TValueType, VRows, VColumns> & inputMatrix)
{
  return inputMatrix.data_block();
};
template <typename TValueType>
const TValueType *
GetPointerToMatrixData(const vnl_matrix<TValueType> & inputMatrix)
{
  return inputMatrix.data_block();
};

template <typename TValueType, unsigned int VRows, unsigned int VColumns>
const TValueType *
GetPointerToMatrixData(const itk::Matrix<TValueType, VRows, VColumns> & inputMatrix)
{
  return inputMatrix.GetVnlMatrix().data_block();
};

/** Sort input to be ordered by magnitude, and returns container with the
 * permutations required for the sorting.
 *
 * For example, if input eigenValues = {10, 0, 40}, the output would be: {2,0,1}
 * and the eigenValues would be modified in-place: {40, 10, 0}.
 *
 * The permutations indices is used to order the matrix of eigenVectors.
 * \sa permuteEigenVectorsWithSortPermutations
 *
 * @tparam TArray  array type with operator []
 * @param eigenValues input array, requires operator []
 * @param numberOfElements size of array
 *
 * @return the permutations needed to sort the input array
 */
template <typename TArray>
std::vector<int>
sortEigenValuesByMagnitude(TArray & eigenValues, const unsigned int numberOfElements)
{
  std::vector<int> indicesSortPermutations(numberOfElements, 0);
  std::iota(std::begin(indicesSortPermutations), std::end(indicesSortPermutations), 0);

  std::sort(std::begin(indicesSortPermutations),
            std::end(indicesSortPermutations),
            [&eigenValues](unsigned int a, unsigned int b) {
              return itk::Math::abs(eigenValues[a]) < itk::Math::abs(eigenValues[b]);
            });
  auto tmpCopy = eigenValues;
  for (unsigned int i = 0; i < numberOfElements; ++i)
  {
    eigenValues[i] = tmpCopy[indicesSortPermutations[i]];
  }
  return indicesSortPermutations;
}

/** Permute a eigenVectors matrix according to the permutation indices
 * computed from the output of a sort function like \sa detail::sortEigenValuesByMagnitude
 *
 * @tparam QMatrix a Eigen3 matrix
 * @param eigenVectors stored in columns
 * @param indicesSortPermutations container with the permutations from the output of
 * a sort function.
 */
template <typename QMatrix>
void
permuteColumnsWithSortIndices(QMatrix & eigenVectors, const std::vector<int> & indicesSortPermutations)
{
  using EigenLibPermutationMatrix = Eigen::PermutationMatrix<Eigen::Dynamic, Eigen::Dynamic>;
  auto numberOfElements = indicesSortPermutations.size();
  // Creates a NxN permutation matrix copying our permutation to the matrix indices.
  // Which holds the 1D array representation of a permutation.
  EigenLibPermutationMatrix perm(numberOfElements);
  perm.setIdentity();
  std::copy(indicesSortPermutations.begin(), indicesSortPermutations.end(), perm.indices().data());
  // Apply it
  eigenVectors = eigenVectors * perm;
}
} // end namespace detail

/** \class SymmetricEigenAnalysisEnums
 * \brief This class contains all enum classes used by SymmetricEigenAnalysis class.
 * \ingroup ITKCommon
 */
class SymmetricEigenAnalysisEnums
{
public:
  /** \class EigenValueOrder
   * \ingroup ITKCommon
   * Order of eigen values
   * OrderByValue:      lambda_1 < lambda_2 < ....
   * OrderByMagnitude:  |lambda_1| < |lambda_2| < .....
   * DoNotOrder:        Default order of eigen values obtained after QL method
   */
  enum class EigenValueOrder : uint8_t
  {
    OrderByValue = 1,
    OrderByMagnitude = 2,
    DoNotOrder = 3
  };
};
// Define how to print enumeration
extern ITKCommon_EXPORT std::ostream &
                        operator<<(std::ostream & out, const SymmetricEigenAnalysisEnums::EigenValueOrder value);

using EigenValueOrderEnum = SymmetricEigenAnalysisEnums::EigenValueOrder;

inline EigenValueOrderEnum
Int2EigenValueOrderEnum(const uint8_t value)
{
  switch (value)
  {
    case 1:
      return EigenValueOrderEnum::OrderByValue;
    case 2:
      return EigenValueOrderEnum::OrderByMagnitude;
    case 3:
      return EigenValueOrderEnum::DoNotOrder;
    default:
      break;
  }
  itkGenericExceptionMacro("Error: Invalid value for conversion.");
}

#if !defined(ITK_LEGACY_REMOVE)
/** Enables reverse compatibility for enumeration values */
static constexpr EigenValueOrderEnum OrderByValue = EigenValueOrderEnum::OrderByValue;
static constexpr EigenValueOrderEnum OrderByMagnitude = EigenValueOrderEnum::OrderByMagnitude;
static constexpr EigenValueOrderEnum DoNotOrder = EigenValueOrderEnum::DoNotOrder;
#endif

/** \class SymmetricEigenAnalysis
 * \brief Find Eigen values of a real 2D symmetric matrix. It
 * serves as a thread-safe alternative to the class:
 * vnl_symmetric_eigensystem, which uses netlib routines.
 *
 * The class is templated over the input matrix (which is expected to provide
 * access to its elements with the [][] operator), matrix to store eigen
 * values (must provide write operations on its elements with the [] operator), and
 * EigenMatrix to store eigen vectors (must provide write access to its elements
 * with the [][] operator).
 *
 * The SetOrderEigenValues() method can be used to order eigen values (and their
 * corresponding eigen vectors if computed) in ascending order. This is the
 * default ordering scheme. Eigen vectors and values can be obtained without
 * ordering by calling SetOrderEigenValues(false).
 *
 * The SetOrderEigenMagnitudes() method can be used to order eigen values (and
 * their corresponding eigen vectors if computed) by magnitude in ascending order.
 *
 * The user of this class is explicitly supposed to set the dimension of the
 * 2D matrix using the SetDimension() method.
 *
 * The class contains routines taken from netlib sources (www.netlib.org).
 * netlib/tql1.c
 * netlib/tql2.c
 * netlib/tred1.c
 * netlib/tred2.c
 *
 * Reference:
 *     num. math. 11, 293-306(1968) by bowdler, martin, reinsch, and
 *     wilkinson.
 *     handbook for auto. comp., vol.ii-linear algebra, 227-240(1971).
 * \ingroup ITKCommon
 */

template <typename TMatrix, typename TVector, typename TEigenMatrix = TMatrix>
class ITK_TEMPLATE_EXPORT SymmetricEigenAnalysis
{
public:
  using EigenValueOrderEnum = itk::EigenValueOrderEnum;
#if !defined(ITK_LEGACY_REMOVE)
  /** Enables reverse compatibility for enumeration values */
  using EigenValueOrderType = EigenValueOrderEnum;
#endif

  SymmetricEigenAnalysis() = default;

  SymmetricEigenAnalysis(const unsigned int dimension)
    : m_Dimension(dimension)
    , m_Order(dimension)
  {}

  ~SymmetricEigenAnalysis() = default;

  using MatrixType = TMatrix;
  using EigenMatrixType = TEigenMatrix;
  using VectorType = TVector;

  /** Compute Eigen values of A
   * A is any type that overloads the [][] operator and contains the
   * symmetric matrix. In practice only the upper triangle of the
   * matrix will be accessed. (Both itk::Matrix and vnl_matrix
   * overload [][] operator.)
   *
   * 'EigenValues' is any type that overloads the [][] operator and will contain
   * the eigen values.
   *
   * No size checking is performed. A is expected to be a square matrix of size
   * m_Dimension.  'EigenValues' is expected to be of length m_Dimension.
   * The matrix is not checked to see if it is symmetric.
   */
  unsigned int
  ComputeEigenValues(const TMatrix & A, TVector & D) const;

  /** Compute Eigen values and vectors of A
   * A is any type that overloads the [][] operator and contains the
   * symmetric matrix. In practice only the upper triangle of the
   * matrix will be accessed. (Both itk::Matrix and vnl_matrix
   * overload [][] operator.)
   *
   * 'EigenValues' is any type that overloads the [] operator and will contain
   * the eigen values.
   *
   * 'EigenVectors' is any type that provides access to its elements with the
   * [][] operator. It is expected be of size m_Dimension * m_Dimension.
   *
   * No size checking is performed. A is expected to be a square matrix of size
   * m_Dimension.  'EigenValues' is expected to be of length m_Dimension.
   * The matrix is not checked to see if it is symmetric.
   *
   * Each row of the matrix 'EigenVectors' represents an eigen vector. (unlike MATLAB
   * where the columns of the [EigenVectors, EigenValues] = eig(A) contains the
   * eigenvectors).
   */
  unsigned int
  ComputeEigenValuesAndVectors(const TMatrix & A, TVector & EigenValues, TEigenMatrix & EigenVectors) const;


  /** Matrix order. Defaults to matrix dimension if not set */
  void
  SetOrder(const unsigned int n)
  {
    m_Order = n;
  }

  /** Get the Matrix order. Will be 0 unless explicitly set, or unless a
   * call to SetDimension has been made in which case it will be the
   * matrix dimension. */
  unsigned int
  GetOrder() const
  {
    return m_Order;
  }

  /** Set/Get methods to order the eigen values in ascending order.
   * This is the default. ie lambda_1 < lambda_2 < ....
   */
  void
  SetOrderEigenValues(const bool b)
  {
    if (b)
    {
      m_OrderEigenValues = EigenValueOrderEnum::OrderByValue;
    }
    else
    {
      m_OrderEigenValues = EigenValueOrderEnum::DoNotOrder;
    }
  }

  bool
  GetOrderEigenValues() const
  {
    return (m_OrderEigenValues == EigenValueOrderEnum::OrderByValue);
  }

  /** Set/Get methods to order the eigen value magnitudes in ascending order.
   * In other words, |lambda_1| < |lambda_2| < .....
   */
  void
  SetOrderEigenMagnitudes(const bool b)
  {
    if (b)
    {
      m_OrderEigenValues = EigenValueOrderEnum::OrderByMagnitude;
    }
    else
    {
      m_OrderEigenValues = EigenValueOrderEnum::DoNotOrder;
    }
  }

  bool
  GetOrderEigenMagnitudes() const
  {
    return (m_OrderEigenValues == EigenValueOrderEnum::OrderByMagnitude);
  }

  /** Set the dimension of the input matrix A. A is a square matrix of
   * size m_Dimension. */
  void
  SetDimension(const unsigned int n)
  {
    m_Dimension = n;
    if (m_Order == 0)
    {
      m_Order = m_Dimension;
    }
  }

  /** Get Matrix dimension, Will be 0 unless explicitly set by a
   * call to SetDimension. */
  unsigned int
  GetDimension() const
  {
    return m_Dimension;
  }

  /** Set/Get to use Eigen library instead of vnl/netlib. */
  void
  SetUseEigenLibrary(const bool input)
  {
    m_UseEigenLibrary = input;
  }
  void
  SetUseEigenLibraryOn()
  {
    m_UseEigenLibrary = true;
  }
  void
  SetUseEigenLibraryOff()
  {
    m_UseEigenLibrary = false;
  }
  bool
  GetUseEigenLibrary() const
  {
    return m_UseEigenLibrary;
  }

private:
  bool                m_UseEigenLibrary{ false };
  unsigned int        m_Dimension{ 0 };
  unsigned int        m_Order{ 0 };
  EigenValueOrderEnum m_OrderEigenValues{ EigenValueOrderEnum::OrderByValue };

  /** Reduces a real symmetric matrix to a symmetric tridiagonal matrix using
   *  orthogonal similarity transformations.
   *  'inputMatrix' contains the real symmetric input matrix. Only the lower
   *  triangle of the matrix need be supplied. The upper triangle is unaltered.
   *  'd' contains the diagonal elements of the tridiagonal matrix.
   *  'e' contains the subdiagonal elements of the tridiagonal matrix in its
   *  last n-1 positions.  e(1) is set to zero.
   *  'e2' contains the squares of the corresponding elements of e.
   *  'e2' may coincide with e if the squares are not needed.
   *  questions and comments should be directed to burton s. garbow.
   *  mathematics and computer science div, argonne national laboratory
   *     this version dated august 1983.
   *
   *  Function adapted from netlib/tred1.c.
   *  [Changed: remove static vars, enforce const correctness.
   *            Use vnl routines as necessary].
   *  Reference:
   *  num. math. 11, 181-195(1968) by martin, reinsch, and wilkinson.
   *    handbook for auto. comp., vol.ii-linear algebra, 212-226(1971).    */
  void
  ReduceToTridiagonalMatrix(double * a, double * d, double * e, double * e2) const;

  /** Reduces a real symmetric matrix to a symmetric tridiagonal matrix using
   *  and accumulating orthogonal similarity transformations.
   *  'inputMatrix' contains the real symmetric input matrix. Only the lower
   *  triangle of the matrix need be supplied. The upper triangle is unaltered.
   *  'diagonalElements' will contains the diagonal elements of the tridiagonal
   *  matrix.
   *  'subDiagonalElements' will contain the subdiagonal elements of the tridiagonal
   *  matrix in its last n-1 positions.  subDiagonalElements(1) is set to zero.
   *  'transformMatrix' contains the orthogonal transformation matrix produced
   *  in the reduction.
   *
   *  Questions and comments should be directed to Burton s. Garbow,
   *  Mathematics and Computer Science Div., Argonne National Laboratory.
   *  This version dated august 1983.
   *
   *  Function adapted from netlib/tred2.c.
   *  [Changed: remove static vars, enforce const correctness.
   *            Use vnl routines as necessary].
   *  Reference:
   *  num. math. 11, 181-195(1968) by martin, reinsch, and wilkinson.
   *    handbook for auto. comp., vol.ii-linear algebra, 212-226(1971).    */
  void
  ReduceToTridiagonalMatrixAndGetTransformation(const double * a, double * d, double * e, double * z) const;

  /** Finds the eigenvalues of a symmetric tridiagonal matrix by the ql method.
   *
   * On input:
   * 'd' contains the diagonal elements of the input matrix.
   * 'e' contains the subdiagonal elements of the input matrix
   * in its last n-1 positions.  e(1) is arbitrary.
   * On Output:
   * 'd' contains the eigenvalues.
   * 'e' has been destroyed.
   *
   * Returns:
   *          zero       for normal return,
   *          j          if the j-th eigenvalue has not been
   *                     determined after 30 iterations.
   *
   *
   * Reference
   *  This subroutine is a translation of the algol procedure tql1,
   *  num. math. 11, 293-306(1968) by bowdler, martin, reinsch, and
   *  wilkinson.
   *  handbook for auto. comp., vol.ii-linear algebra, 227-240(1971).
   *
   *  Questions and comments should be directed to Burton s. Garbow,
   *  Mathematics and Computer Science Div., Argonne National Laboratory.
   *  This version dated august 1983.
   *
   *  Function Adapted from netlib/tql1.c.
   *  [Changed: remove static vars, enforce const correctness.
   *            Use vnl routines as necessary]                      */
  unsigned int
  ComputeEigenValuesUsingQL(double * d, double * e) const;

  /** Finds the eigenvalues and eigenvectors of a symmetric tridiagonal matrix
   * by the ql method.
   *
   * On input:
   * 'd' contains the diagonal elements of the input matrix.
   * 'e' contains the subdiagonal elements of the input matrix
   * in its last n-1 positions.  e(1) is arbitrary.
   * 'z' contains the transformation matrix produced in the reduction by
   * ReduceToTridiagonalMatrixAndGetTransformation(), if performed. If the
   * eigenvectors of the tridiagonal matrix are desired, z must contain
   * the identity matrix.

   * On Output:
   * 'd' contains the eigenvalues.
   * 'e' has been destroyed.
   * 'z' contains orthonormal eigenvectors of the symmetric tridiagonal
   * (or full) matrix.
   *
   * Returns:
   *          zero       for normal return,
   *          j          if the j-th eigenvalue has not been
   *                     determined after 1000 iterations.
   *
   * Reference
   *  This subroutine is a translation of the algol procedure tql1,
   *  num. math. 11, 293-306(1968) by bowdler, martin, reinsch, and
   *  wilkinson.
   *  handbook for auto. comp., vol.ii-linear algebra, 227-240(1971).
   *
   *  Questions and comments should be directed to Burton s. Garbow,
   *  Mathematics and Computer Science Div., Argonne National Laboratory.
   *  This version dated august 1983.
   *
   *  Function Adapted from netlib/tql2.c.
   *  [Changed: remove static vars, enforce const correctness.
   *            Use vnl routines as necessary]
   */
  unsigned int
  ComputeEigenValuesAndVectorsUsingQL(double * d, double * e, double * z) const;

  /* Legacy algorithms using thread-safe netlib.
   * \sa ComputeEigenValues and \sa ComputeEigenValuesAndVectors
   */
  unsigned int
  ComputeEigenValuesLegacy(const TMatrix & A, TVector & D) const;

  unsigned int
  ComputeEigenValuesAndVectorsLegacy(const TMatrix & A, TVector & EigenValues, TEigenMatrix & EigenVectors) const;

  /* Helper to get the matrix value type for EigenLibMatrix typename.
   *
   * If the TMatrix is vnl, the type is in element_type.
   * In TMatrix is itk::Matrix, or any itk::FixedArray is in ValueType.
   *
   * To use this function:
   * using ValueType = decltype(this->GetMatrixType(true));
   *
   * \note The two `GetMatrixValueType` overloads have different
   * parameter declarations (`bool` and `...`), to avoid that both
   * functions are equally good candidates during overload resolution,
   * in case `element_type` and `ValueType` are both nested types of
   * `TMatrix` (which is the case when `TMatrix` = `itk::Array2D`).
   */
  template <typename QMatrix = TMatrix>
  auto
  GetMatrixValueType(bool) const -> typename QMatrix::element_type
  {
    return QMatrix::element_type();
  }
  template <typename QMatrix = TMatrix>
  auto
  GetMatrixValueType(...) const -> typename QMatrix::ValueType
  {
    return QMatrix::ValueType();
  }

  /* Wrapper that call the right implementation for the type of matrix.  */
  unsigned int
  ComputeEigenValuesAndVectorsWithEigenLibrary(const TMatrix & A,
                                               TVector &       EigenValues,
                                               TEigenMatrix &  EigenVectors) const
  {
    return ComputeEigenValuesAndVectorsWithEigenLibraryImpl(A, EigenValues, EigenVectors, true);
  }

  /* Implementation detail using EigenLib that performs a copy of the input matrix.
   *
   * @param (long) implementation detail argument making this implementation less favourable
   *   to be chosen if alternatives are available.
   *
   * @return an unsigned int with no information value (no error code in EigenLib) */
  template <typename QMatrix>
  auto
  ComputeEigenValuesAndVectorsWithEigenLibraryImpl(const QMatrix & A,
                                                   TVector &       EigenValues,
                                                   TEigenMatrix &  EigenVectors,
                                                   long) const -> decltype(1U)
  {
    using ValueType = decltype(GetMatrixValueType(true));
    using EigenLibMatrixType = Eigen::Matrix<ValueType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>;
    EigenLibMatrixType inputMatrix(m_Dimension, m_Dimension);
    for (unsigned int row = 0; row < m_Dimension; ++row)
    {
      for (unsigned int col = 0; col < m_Dimension; ++col)
      {
        inputMatrix(row, col) = A(row, col);
      }
    }
    using EigenSolverType = Eigen::SelfAdjointEigenSolver<EigenLibMatrixType>;
    EigenSolverType solver(inputMatrix); // Computes EigenValues and EigenVectors
    const auto &    eigenValues = solver.eigenvalues();
    /* Column  k  of the returned matrix is an eigenvector corresponding to
     * eigenvalue number $ k $ as returned by eigenvalues().
     * The eigenvectors are normalized to have (Euclidean) norm equal to one. */
    const auto & eigenVectors = solver.eigenvectors();

    if (m_OrderEigenValues == EigenValueOrderEnum::OrderByMagnitude)
    {
      auto copyEigenValues = eigenValues;
      auto copyEigenVectors = eigenVectors;
      auto indicesSortPermutations = detail::sortEigenValuesByMagnitude(copyEigenValues, m_Dimension);
      detail::permuteColumnsWithSortIndices(copyEigenVectors, indicesSortPermutations);

      for (unsigned int row = 0; row < m_Dimension; ++row)
      {
        EigenValues[row] = copyEigenValues[row];
        for (unsigned int col = 0; col < m_Dimension; ++col)
        {
          EigenVectors[row][col] = copyEigenVectors(col, row);
        }
      }
    }
    else
    {
      for (unsigned int row = 0; row < m_Dimension; ++row)
      {
        EigenValues[row] = eigenValues[row];
        for (unsigned int col = 0; col < m_Dimension; ++col)
        {
          EigenVectors[row][col] = eigenVectors(col, row);
        }
      }
    }
    // No error code
    return 1;
  }


  /* Implementation detail using EigenLib that do not perform a copy.
   * It needs the existence of a pointer to matrix data. \sa GetPointerToMatrixData
   * If new types want to use this method, an appropriate overload of GetPointerToMatrixData
   * should be included.
   *
   * @param (bool) implementation detail argument making this implementation the most favourable
   *   to be chosen from all the alternative implementations.
   *
   * @return an unsigned int with no information value (no error code in EigenLib) */
  template <typename QMatrix>
  auto
  ComputeEigenValuesAndVectorsWithEigenLibraryImpl(const QMatrix & A,
                                                   TVector &       EigenValues,
                                                   TEigenMatrix &  EigenVectors,
                                                   bool) const -> decltype(GetPointerToMatrixData(A), 1U)
  {
    auto pointerToData = GetPointerToMatrixData(A);
    using PointerType = decltype(pointerToData);
    using ValueTypeCV = std::remove_pointer_t<PointerType>;
    using ValueType = std::remove_cv_t<ValueTypeCV>;
    using EigenLibMatrixType = Eigen::Matrix<ValueType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>;
    using EigenConstMatrixMap = Eigen::Map<const EigenLibMatrixType>;
    EigenConstMatrixMap inputMatrix(pointerToData, m_Dimension, m_Dimension);
    using EigenSolverType = Eigen::SelfAdjointEigenSolver<EigenLibMatrixType>;
    EigenSolverType solver(inputMatrix); // Computes EigenValues and EigenVectors
    const auto &    eigenValues = solver.eigenvalues();
    /* Column  k  of the returned matrix is an eigenvector corresponding to
     * eigenvalue number $ k $ as returned by eigenvalues().
     * The eigenvectors are normalized to have (Euclidean) norm equal to one. */
    const auto & eigenVectors = solver.eigenvectors();
    if (m_OrderEigenValues == EigenValueOrderEnum::OrderByMagnitude)
    {
      auto copyEigenValues = eigenValues;
      auto copyEigenVectors = eigenVectors;
      auto indicesSortPermutations = detail::sortEigenValuesByMagnitude(copyEigenValues, m_Dimension);
      detail::permuteColumnsWithSortIndices(copyEigenVectors, indicesSortPermutations);
      for (unsigned int row = 0; row < m_Dimension; ++row)
      {
        EigenValues[row] = copyEigenValues[row];
        for (unsigned int col = 0; col < m_Dimension; ++col)
        {
          EigenVectors[row][col] = copyEigenVectors(col, row);
        }
      }
    }
    else
    {
      for (unsigned int row = 0; row < m_Dimension; ++row)
      {
        EigenValues[row] = eigenValues[row];
        for (unsigned int col = 0; col < m_Dimension; ++col)
        {
          EigenVectors[row][col] = eigenVectors(col, row);
        }
      }
    }
    // No error code
    return 1;
  }

  /* Wrapper that call the right implementation for the type of matrix.  */
  unsigned int
  ComputeEigenValuesWithEigenLibrary(const TMatrix & A, TVector & EigenValues) const
  {
    return ComputeEigenValuesWithEigenLibraryImpl(A, EigenValues, true);
  }

  /* Implementation detail using EigenLib that performs a copy of the input matrix.
   *
   * @param (long) implementation detail argument making this implementation less favourable
   *   to be chosen if alternatives are available.
   *
   * @return an unsigned int with no information value (no error code in EigenLib) */
  template <typename QMatrix>
  auto
  ComputeEigenValuesWithEigenLibraryImpl(const QMatrix & A, TVector & EigenValues, long) const -> decltype(1U)
  {
    using ValueType = decltype(GetMatrixValueType(true));
    using EigenLibMatrixType = Eigen::Matrix<ValueType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>;
    EigenLibMatrixType inputMatrix(m_Dimension, m_Dimension);
    for (unsigned int row = 0; row < m_Dimension; ++row)
    {
      for (unsigned int col = 0; col < m_Dimension; ++col)
      {
        inputMatrix(row, col) = A(row, col);
      }
    }
    using EigenSolverType = Eigen::SelfAdjointEigenSolver<EigenLibMatrixType>;
    EigenSolverType solver(inputMatrix, Eigen::EigenvaluesOnly);
    auto            eigenValues = solver.eigenvalues();
    if (m_OrderEigenValues == EigenValueOrderEnum::OrderByMagnitude)
    {
      detail::sortEigenValuesByMagnitude(eigenValues, m_Dimension);
    }
    for (unsigned int i = 0; i < m_Dimension; ++i)
    {
      EigenValues[i] = eigenValues[i];
    }

    // No error code
    return 1;
  }

  /* Implementation detail using EigenLib that do not perform a copy.
   * It needs the existence of a pointer to matrix data. \sa GetPointerToMatrixData
   * If new types want to use this method, an appropriate overload of GetPointerToMatrixData
   * should be included.
   *
   * @param (bool) implementation detail argument making this implementation the most favourable
   *   to be chosen from all the alternative implementations.
   *
   * @return an unsigned int with no information value (no error code in EigenLib) */
  template <typename QMatrix>
  auto
  ComputeEigenValuesWithEigenLibraryImpl(const QMatrix & A, TVector & EigenValues, bool) const
    -> decltype(GetPointerToMatrixData(A), 1U)
  {
    auto pointerToData = GetPointerToMatrixData(A);
    using PointerType = decltype(pointerToData);
    using ValueTypeCV = std::remove_pointer_t<PointerType>;
    using ValueType = std::remove_cv_t<ValueTypeCV>;
    using EigenLibMatrixType = Eigen::Matrix<ValueType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>;
    using EigenConstMatrixMap = Eigen::Map<const EigenLibMatrixType>;
    EigenConstMatrixMap inputMatrix(pointerToData, m_Dimension, m_Dimension);
    using EigenSolverType = Eigen::SelfAdjointEigenSolver<EigenLibMatrixType>;
    EigenSolverType solver(inputMatrix, Eigen::EigenvaluesOnly);
    auto            eigenValues = solver.eigenvalues();
    if (m_OrderEigenValues == EigenValueOrderEnum::OrderByMagnitude)
    {
      detail::sortEigenValuesByMagnitude(eigenValues, m_Dimension);
    }
    for (unsigned int i = 0; i < m_Dimension; ++i)
    {
      EigenValues[i] = eigenValues[i];
    }
    // No error code
    return 1;
  }
};

template <typename TMatrix, typename TVector, typename TEigenMatrix>
std::ostream &
operator<<(std::ostream & os, const SymmetricEigenAnalysis<TMatrix, TVector, TEigenMatrix> & s)
{
  os << "[ClassType: SymmetricEigenAnalysis]" << std::endl;
  os << "  Dimension : " << s.GetDimension() << std::endl;
  os << "  Order : " << s.GetOrder() << std::endl;
  os << "  OrderEigenValues: " << s.GetOrderEigenValues() << std::endl;
  os << "  OrderEigenMagnitudes: " << s.GetOrderEigenMagnitudes() << std::endl;
  os << "  UseEigenLibrary: " << s.GetUseEigenLibrary() << std::endl;
  return os;
}

template <unsigned int VDimension, typename TMatrix, typename TVector, typename TEigenMatrix = TMatrix>
class ITK_TEMPLATE_EXPORT SymmetricEigenAnalysisFixedDimension
{
public:
#if !defined(ITK_LEGACY_REMOVE)
  /** Enables reverse compatibility for enumeration values */
  using EigenValueOrderType = EigenValueOrderEnum;
#endif
#if !defined(ITK_LEGACY_REMOVE)
  // We need to expose the enum values at the class level
  // for backwards compatibility
  static constexpr EigenValueOrderEnum OrderByValue = EigenValueOrderEnum::OrderByValue;
  static constexpr EigenValueOrderEnum OrderByMagnitude = EigenValueOrderEnum::OrderByMagnitude;
  static constexpr EigenValueOrderEnum DoNotOrder = EigenValueOrderEnum::DoNotOrder;
#endif

  using MatrixType = TMatrix;
  using EigenMatrixType = TEigenMatrix;
  using VectorType = TVector;

  /** Compute Eigen values of A
   * A is any type that overloads the [][] operator and contains the
   * symmetric matrix. In practice only the upper triangle of the
   * matrix will be accessed. (Both itk::Matrix and vnl_matrix
   * overload [][] operator.)
   *
   * 'EigenValues' is any type that overloads the [] operator and will contain
   * the eigen values.
   *
   * No size checking is performed. A is expected to be a square matrix of size
   * VDimension.  'EigenValues' is expected to be of length VDimension.
   * The matrix is not checked to see if it is symmetric.
   */
  unsigned int
  ComputeEigenValues(const TMatrix & A, TVector & EigenValues) const
  {
    return ComputeEigenValuesWithEigenLibraryImpl(A, EigenValues, true);
  }

  /** Compute Eigen values and vectors of A
   * A is any type that overloads the [][] operator and contains the
   * symmetric matrix. In practice only the upper triangle of the
   * matrix will be accessed. (Both itk::Matrix and vnl_matrix
   * overload [][] operator.)
   *
   * 'EigenValues' is any type that overloads the [] operator and will contain
   * the eigen values.
   *
   * 'EigenVectors' is any type that provides access to its elements with the
   * [][] operator. It is expected be of size VDimension * VDimension.
   *
   * No size checking is performed. A is expected to be a square matrix of size
   * VDimension.  'EigenValues' is expected to be of length VDimension.
   * The matrix is not checked to see if it is symmetric.
   *
   * Each row of the matrix 'EigenVectors' represents an eigen vector. (unlike MATLAB
   * where the columns of the [EigenVectors, EigenValues] = eig(A) contains the
   * eigenvectors).
   */
  unsigned int
  ComputeEigenValuesAndVectors(const TMatrix & A, TVector & EigenValues, TEigenMatrix & EigenVectors) const
  {
    return ComputeEigenValuesAndVectorsWithEigenLibraryImpl(A, EigenValues, EigenVectors, true);
  }

  void
  SetOrderEigenValues(const bool b)
  {
    if (b)
    {
      m_OrderEigenValues = EigenValueOrderEnum::OrderByValue;
    }
    else
    {
      m_OrderEigenValues = EigenValueOrderEnum::DoNotOrder;
    }
  }
  bool
  GetOrderEigenValues() const
  {
    return (m_OrderEigenValues == EigenValueOrderEnum::OrderByValue);
  }
  void
  SetOrderEigenMagnitudes(const bool b)
  {
    if (b)
    {
      m_OrderEigenValues = EigenValueOrderEnum::OrderByMagnitude;
    }
    else
    {
      m_OrderEigenValues = EigenValueOrderEnum::DoNotOrder;
    }
  }
  bool
  GetOrderEigenMagnitudes() const
  {
    return (m_OrderEigenValues == EigenValueOrderEnum::OrderByMagnitude);
  }
  constexpr unsigned int
  GetOrder() const
  {
    return VDimension;
  }
  constexpr unsigned int
  GetDimension() const
  {
    return VDimension;
  }
  constexpr bool
  GetUseEigenLibrary() const
  {
    return true;
  }

private:
  EigenValueOrderEnum m_OrderEigenValues{ EigenValueOrderEnum::OrderByValue };

  /* Helper to get the matrix value type for EigenLibMatrix typename.
   *
   * If the TMatrix is vnl, the type is in element_type.
   * In TMatrix is itk::Matrix, or any itk::FixedArray is in ValueType.
   *
   * To use this function:
   * using ValueType = decltype(this->GetMatrixType(true));
   */
  template <typename QMatrix = TMatrix>
  auto
  GetMatrixValueType(bool) const -> typename QMatrix::element_type
  {
    return QMatrix::element_type();
  }
  template <typename QMatrix = TMatrix>
  auto
  GetMatrixValueType(bool) const -> typename QMatrix::ValueType
  {
    return QMatrix::ValueType();
  }

  /* Implementation detail using EigenLib that do not perform a copy.
   * It needs the existence of a pointer to matrix data. \sa GetPointerToMatrixData
   * If new types want to use this method, an appropriate overload of GetPointerToMatrixData
   * should be included.
   *
   * @param (bool) implementation detail argument making this implementation the most favourable
   *   to be chosen from all the alternative implementations.
   *
   * @return an unsigned int with no information value (no error code in EigenLib) */
  template <typename QMatrix>
  auto
  ComputeEigenValuesAndVectorsWithEigenLibraryImpl(const QMatrix & A,
                                                   TVector &       EigenValues,
                                                   TEigenMatrix &  EigenVectors,
                                                   bool) const -> decltype(GetPointerToMatrixData(A), 1U)
  {
    auto pointerToData = GetPointerToMatrixData(A);
    using PointerType = decltype(pointerToData);
    using ValueTypeCV = std::remove_pointer_t<PointerType>;
    using ValueType = std::remove_cv_t<ValueTypeCV>;
    using EigenLibMatrixType = Eigen::Matrix<ValueType, VDimension, VDimension, Eigen::RowMajor>;
    using EigenConstMatrixMap = Eigen::Map<const EigenLibMatrixType>;
    EigenConstMatrixMap inputMatrix(pointerToData);
    using EigenSolverType = Eigen::SelfAdjointEigenSolver<EigenLibMatrixType>;
    EigenSolverType solver(inputMatrix); // Computes EigenValues and EigenVectors
    const auto &    eigenValues = solver.eigenvalues();
    /* Column  k  of the returned matrix is an eigenvector corresponding to
     * eigenvalue number $ k $ as returned by eigenvalues().
     * The eigenvectors are normalized to have (Euclidean) norm equal to one. */
    const auto & eigenVectors = solver.eigenvectors();
    if (m_OrderEigenValues == EigenValueOrderEnum::OrderByMagnitude)
    {
      auto copyEigenValues = eigenValues;
      auto copyEigenVectors = eigenVectors;
      auto indicesSortPermutations = detail::sortEigenValuesByMagnitude(copyEigenValues, VDimension);
      detail::permuteColumnsWithSortIndices(copyEigenVectors, indicesSortPermutations);
      for (unsigned int row = 0; row < VDimension; ++row)
      {
        EigenValues[row] = copyEigenValues[row];
        for (unsigned int col = 0; col < VDimension; ++col)
        {
          EigenVectors[row][col] = copyEigenVectors(col, row);
        }
      }
    }
    else
    {
      for (unsigned int row = 0; row < VDimension; ++row)
      {
        EigenValues[row] = eigenValues[row];
        for (unsigned int col = 0; col < VDimension; ++col)
        {
          EigenVectors[row][col] = eigenVectors(col, row);
        }
      }
    }
    // No error code
    return 1;
  }

  /* Implementation detail using EigenLib that performs a copy of the input matrix.
   *
   * @param (long) implementation detail argument making this implementation less favourable
   *   to be chosen if alternatives are available.
   *
   * @return an unsigned int with no information value (no error code in EigenLib) */
  template <typename QMatrix>
  auto
  ComputeEigenValuesAndVectorsWithEigenLibraryImpl(const QMatrix & A,
                                                   TVector &       EigenValues,
                                                   TEigenMatrix &  EigenVectors,
                                                   long) const -> decltype(1U)
  {
    using ValueType = decltype(GetMatrixValueType(true));
    using EigenLibMatrixType = Eigen::Matrix<ValueType, VDimension, VDimension, Eigen::RowMajor>;
    EigenLibMatrixType inputMatrix;
    for (unsigned int row = 0; row < VDimension; ++row)
    {
      for (unsigned int col = 0; col < VDimension; ++col)
      {
        inputMatrix(row, col) = A(row, col);
      }
    }
    using EigenSolverType = Eigen::SelfAdjointEigenSolver<EigenLibMatrixType>;
    EigenSolverType solver(inputMatrix); // Computes EigenValues and EigenVectors
    const auto &    eigenValues = solver.eigenvalues();
    /* Column  k  of the returned matrix is an eigenvector corresponding to
     * eigenvalue number $ k $ as returned by eigenvalues().
     * The eigenvectors are normalized to have (Euclidean) norm equal to one. */
    const auto & eigenVectors = solver.eigenvectors();

    if (m_OrderEigenValues == EigenValueOrderEnum::OrderByMagnitude)
    {
      auto copyEigenValues = eigenValues;
      auto copyEigenVectors = eigenVectors;
      auto indicesSortPermutations = detail::sortEigenValuesByMagnitude(copyEigenValues, VDimension);
      detail::permuteColumnsWithSortIndices(copyEigenVectors, indicesSortPermutations);

      for (unsigned int row = 0; row < VDimension; ++row)
      {
        EigenValues[row] = copyEigenValues[row];
        for (unsigned int col = 0; col < VDimension; ++col)
        {
          EigenVectors[row][col] = copyEigenVectors(col, row);
        }
      }
    }
    else
    {
      for (unsigned int row = 0; row < VDimension; ++row)
      {
        EigenValues[row] = eigenValues[row];
        for (unsigned int col = 0; col < VDimension; ++col)
        {
          EigenVectors[row][col] = eigenVectors(col, row);
        }
      }
    }
    // No error code
    return 1;
  }

  /* Implementation detail using EigenLib that performs a copy of the input matrix.
   *
   * @param (long) implementation detail argument making this implementation less favourable
   *   to be chosen if alternatives are available.
   *
   * @return an unsigned int with no information value (no error code in EigenLib) */
  template <typename QMatrix>
  auto
  ComputeEigenValuesWithEigenLibraryImpl(const QMatrix & A, TVector & EigenValues, long) const -> decltype(1U)
  {
    using ValueType = decltype(GetMatrixValueType(true));
    using EigenLibMatrixType = Eigen::Matrix<ValueType, VDimension, VDimension, Eigen::RowMajor>;
    EigenLibMatrixType inputMatrix;
    for (unsigned int row = 0; row < VDimension; ++row)
    {
      for (unsigned int col = 0; col < VDimension; ++col)
      {
        inputMatrix(row, col) = A(row, col);
      }
    }
    using EigenSolverType = Eigen::SelfAdjointEigenSolver<EigenLibMatrixType>;
    EigenSolverType solver(inputMatrix, Eigen::EigenvaluesOnly);
    auto            eigenValues = solver.eigenvalues();
    if (m_OrderEigenValues == EigenValueOrderEnum::OrderByMagnitude)
    {
      detail::sortEigenValuesByMagnitude(eigenValues, VDimension);
    }
    for (unsigned int i = 0; i < VDimension; ++i)
    {
      EigenValues[i] = eigenValues[i];
    }

    // No error code
    return 1;
  }

  /* Implementation detail using EigenLib that do not perform a copy.
   * It needs the existence of a pointer to matrix data. \sa GetPointerToMatrixData
   * If new types want to use this method, an appropriate overload of GetPointerToMatrixData
   * should be included.
   *
   * @param (bool) implementation detail argument making this implementation the most favourable
   *   to be chosen from all the alternative implementations.
   *
   * @return an unsigned int with no information value (no error code in EigenLib) */
  template <typename QMatrix>
  auto
  ComputeEigenValuesWithEigenLibraryImpl(const QMatrix & A, TVector & EigenValues, bool) const
    -> decltype(GetPointerToMatrixData(A), 1U)
  {
    auto pointerToData = GetPointerToMatrixData(A);
    using PointerType = decltype(pointerToData);
    using ValueTypeCV = std::remove_pointer_t<PointerType>;
    using ValueType = std::remove_cv_t<ValueTypeCV>;
    using EigenLibMatrixType = Eigen::Matrix<ValueType, VDimension, VDimension, Eigen::RowMajor>;
    using EigenConstMatrixMap = Eigen::Map<const EigenLibMatrixType>;
    EigenConstMatrixMap inputMatrix(pointerToData);
    using EigenSolverType = Eigen::SelfAdjointEigenSolver<EigenLibMatrixType>;
    EigenSolverType solver(inputMatrix, Eigen::EigenvaluesOnly);
    auto            eigenValues = solver.eigenvalues();
    if (m_OrderEigenValues == EigenValueOrderEnum::OrderByMagnitude)
    {
      detail::sortEigenValuesByMagnitude(eigenValues, VDimension);
    }
    for (unsigned int i = 0; i < VDimension; ++i)
    {
      EigenValues[i] = eigenValues[i];
    }
    // No error code
    return 1;
  }
};

template <unsigned int VDimension, typename TMatrix, typename TVector, typename TEigenMatrix>
std::ostream &
operator<<(std::ostream &                                                                           os,
           const SymmetricEigenAnalysisFixedDimension<VDimension, TMatrix, TVector, TEigenMatrix> & s)
{
  os << "[ClassType: SymmetricEigenAnalysisFixedDimension]" << std::endl;
  os << "  Dimension : " << s.GetDimension() << std::endl;
  os << "  Order : " << s.GetOrder() << std::endl;
  os << "  OrderEigenValues: " << s.GetOrderEigenValues() << std::endl;
  os << "  OrderEigenMagnitudes: " << s.GetOrderEigenMagnitudes() << std::endl;
  os << "  UseEigenLibrary: " << s.GetUseEigenLibrary() << std::endl;
  return os;
}
} // end namespace itk

#ifndef ITK_MANUAL_INSTANTIATION
#  include "itkSymmetricEigenAnalysis.hxx"
#endif

#endif