Blaze 3.9
gesdd.h
Go to the documentation of this file.
1//=================================================================================================
33//=================================================================================================
34
35#ifndef _BLAZE_MATH_LAPACK_GESDD_H_
36#define _BLAZE_MATH_LAPACK_GESDD_H_
37
38
39//*************************************************************************************************
40// Includes
41//*************************************************************************************************
42
43#include <memory>
44#include <blaze/math/Aliases.h>
57#include <blaze/util/Assert.h>
58#include <blaze/util/EnableIf.h>
61
62
63namespace blaze {
64
65//=================================================================================================
66//
67// LAPACK SVD FUNCTIONS (GESDD)
68//
69//=================================================================================================
70
71//*************************************************************************************************
74template< typename MT, bool SO, typename VT, bool TF >
75void gesdd( DenseMatrix<MT,SO>& A, DenseVector<VT,TF>& s );
76
77template< typename MT1, bool SO, typename MT2, typename VT, bool TF >
78void gesdd( DenseMatrix<MT1,SO>& A, DenseMatrix<MT2,SO>& U, DenseVector<VT,TF>& s, char jobz );
79
80template< typename MT1, bool SO, typename MT2, typename VT, bool TF >
81void gesdd( DenseMatrix<MT1,SO>& A, DenseVector<VT,TF>& s, DenseMatrix<MT2,SO>& V, char jobz );
82
83template< typename MT1, bool SO, typename MT2, typename VT, bool TF, typename MT3 >
84void gesdd( DenseMatrix<MT1,SO>& A, DenseMatrix<MT2,SO>& U,
85 DenseVector<VT,TF>& s, DenseMatrix<MT3,SO>& V, char jobz );
87//*************************************************************************************************
88
89
90//*************************************************************************************************
107template< typename MT // Type of the matrix A
108 , bool SO // Storage order of the matrix A
109 , typename VT // Type of the vector s
110 , bool TF > // Transpose flag of the vector s
111inline auto gesdd_backend( DenseMatrix<MT,SO>& A, DenseVector<VT,TF>& s )
112 -> DisableIf_t< IsComplex_v< ElementType_t<MT> > >
113{
114 BLAZE_INTERNAL_ASSERT( (*s).size() == min( (*A).rows(), (*A).columns() ), "Invalid vector dimension detected" );
115
116 using ET = ElementType_t<MT>;
117
118 blas_int_t m ( numeric_cast<blas_int_t>( SO ? (*A).rows() : (*A).columns() ) );
119 blas_int_t n ( numeric_cast<blas_int_t>( SO ? (*A).columns() : (*A).rows() ) );
120 blas_int_t lda ( numeric_cast<blas_int_t>( (*A).spacing() ) );
121 blas_int_t info( 0 );
122
123 const blas_int_t minimum( min( m, n ) );
124 const blas_int_t maximum( max( m, n ) );
125
126 blas_int_t lwork( 3*minimum + max( maximum, 7*minimum ) + 2 );
127 const std::unique_ptr<ET[]> work( new ET[lwork] );
128 const std::unique_ptr<blas_int_t[]> iwork( new blas_int_t[8*minimum] );
129
130 gesdd( 'N', m, n, (*A).data(), lda, (*s).data(),
131 nullptr, 1, nullptr, 1, work.get(), lwork, iwork.get(), &info );
132
133 BLAZE_INTERNAL_ASSERT( info >= 0, "Invalid argument for singular value decomposition" );
134
135 if( info > 0 ) {
136 BLAZE_THROW_LAPACK_ERROR( "Singular value decomposition failed" );
137 }
138}
140//*************************************************************************************************
141
142
143//*************************************************************************************************
160template< typename MT // Type of the matrix A
161 , bool SO // Storage order of the matrix A
162 , typename VT // Type of the vector s
163 , bool TF > // Transpose flag of the vector s
164inline auto gesdd_backend( DenseMatrix<MT,SO>& A, DenseVector<VT,TF>& s )
165 -> EnableIf_t< IsComplex_v< ElementType_t<MT> > >
166{
167 BLAZE_INTERNAL_ASSERT( (*s).size() == min( (*A).rows(), (*A).columns() ), "Invalid vector dimension detected" );
168
169 using CT = ElementType_t<MT>;
170 using BT = UnderlyingElement_t<CT>;
171
172 blas_int_t m ( numeric_cast<blas_int_t>( SO ? (*A).rows() : (*A).columns() ) );
173 blas_int_t n ( numeric_cast<blas_int_t>( SO ? (*A).columns() : (*A).rows() ) );
174 blas_int_t lda ( numeric_cast<blas_int_t>( (*A).spacing() ) );
175 blas_int_t info( 0 );
176
177 const blas_int_t minimum( min( m, n ) );
178 const blas_int_t maximum( max( m, n ) );
179
180 blas_int_t lwork( 2*minimum + maximum + 2 );
181 const std::unique_ptr<CT[]> work( new CT[lwork] );
182 const std::unique_ptr<BT[]> rwork( new BT[7*minimum] );
183 const std::unique_ptr<blas_int_t[]> iwork( new blas_int_t[8*minimum] );
184
185 gesdd( 'N', m, n, (*A).data(), lda, (*s).data(),
186 nullptr, 1, nullptr, 1, work.get(), lwork, rwork.get(), iwork.get(), &info );
187
188 BLAZE_INTERNAL_ASSERT( info >= 0, "Invalid argument for singular value decomposition" );
189
190 if( info > 0 ) {
191 BLAZE_THROW_LAPACK_ERROR( "Singular value decomposition failed" );
192 }
193}
195//*************************************************************************************************
196
197
198//*************************************************************************************************
260template< typename MT // Type of the matrix A
261 , bool SO // Storage order of the matrix A
262 , typename VT // Type of the vector s
263 , bool TF > // Transpose flag of the vector s
265{
271
276
277 const size_t M( (*A).rows() );
278 const size_t N( (*A).columns() );
279 const size_t mindim( min( M, N ) );
280
281 resize( *s, mindim, false );
282
283 if( M == 0UL || N == 0 ) {
284 return;
285 }
286
287 gesdd_backend( A, s );
288}
289//*************************************************************************************************
290
291
292//*************************************************************************************************
311template< typename MT1 // Type of the matrix A
312 , bool SO // Storage order of all matrices
313 , typename MT2 // Type of the matrix U
314 , typename VT // Type of the vector s
315 , bool TF > // Transpose flag of the vector s
316inline auto gesdd_backend( DenseMatrix<MT1,SO>& A, DenseMatrix<MT2,SO>& U,
317 DenseVector<VT,TF>& s, char jobz )
318 -> DisableIf_t< IsComplex_v< ElementType_t<MT1> > >
319{
320 BLAZE_INTERNAL_ASSERT( jobz == 'O' || jobz == 'N', "Invalid jobz flag detected" );
321 BLAZE_INTERNAL_ASSERT( jobz == 'N' || isSquare( *U ), "Invalid matrix dimensions detected" );
322 BLAZE_INTERNAL_ASSERT( jobz == 'N' || (*U).rows() == (*A).rows(), "Invalid matrix dimensions detected" );
323 BLAZE_INTERNAL_ASSERT( (*s).size() == min( (*A).rows(), (*A).columns() ), "Invalid vector dimension detected" );
324
325 using ET = ElementType_t<MT1>;
326
327 blas_int_t m ( numeric_cast<blas_int_t>( SO ? (*A).rows() : (*A).columns() ) );
328 blas_int_t n ( numeric_cast<blas_int_t>( SO ? (*A).columns() : (*A).rows() ) );
329 blas_int_t lda ( numeric_cast<blas_int_t>( (*A).spacing() ) );
330 blas_int_t ldu ( numeric_cast<blas_int_t>( (*U).spacing() ) );
331 blas_int_t info( 0 );
332
333 const blas_int_t minimum( min( m, n ) );
334 const blas_int_t maximum( max( m, n ) );
335
336 blas_int_t lwork( ( jobz == 'O' )
337 ?( 3*minimum + max( maximum, 5*minimum*minimum + 4*maximum ) + 2 )
338 :( 3*minimum + max( maximum, 7*minimum ) + 2 ) );
339 const std::unique_ptr<ET[]> work( new ET[lwork] );
340 const std::unique_ptr<blas_int_t[]> iwork( new blas_int_t[8*minimum] );
341
342 gesdd( jobz, m, n, (*A).data(), lda, (*s).data(),
343 ( SO ? (*U).data() : nullptr ), ( SO ? ldu : 1 ),
344 ( SO ? nullptr : (*U).data() ), ( SO ? 1 : ldu ),
345 work.get(), lwork, iwork.get(), &info );
346
347 BLAZE_INTERNAL_ASSERT( info >= 0, "Invalid argument for singular value decomposition" );
348
349 if( info > 0 ) {
350 BLAZE_THROW_LAPACK_ERROR( "Singular value decomposition failed" );
351 }
352}
354//*************************************************************************************************
355
356
357//*************************************************************************************************
376template< typename MT1 // Type of the matrix A
377 , bool SO // Storage order of all matrices
378 , typename MT2 // Type of the matrix U
379 , typename VT // Type of the vector s
380 , bool TF > // Transpose flag of the vector s
381inline auto gesdd_backend( DenseMatrix<MT1,SO>& A, DenseMatrix<MT2,SO>& U,
382 DenseVector<VT,TF>& s, char jobz )
383 -> EnableIf_t< IsComplex_v< ElementType_t<MT1> > >
384{
385 BLAZE_INTERNAL_ASSERT( jobz == 'O' || jobz == 'N', "Invalid jobz flag detected" );
386 BLAZE_INTERNAL_ASSERT( jobz == 'N' || isSquare( *U ), "Invalid matrix dimensions detected" );
387 BLAZE_INTERNAL_ASSERT( jobz == 'N' || (*U).rows() == (*A).rows(), "Invalid matrix dimensions detected" );
388 BLAZE_INTERNAL_ASSERT( (*s).size() == min( (*A).rows(), (*A).columns() ), "Invalid vector dimension detected" );
389
390 using CT = ElementType_t<MT1>;
391 using BT = UnderlyingElement_t<CT>;
392
393 blas_int_t m ( numeric_cast<blas_int_t>( SO ? (*A).rows() : (*A).columns() ) );
394 blas_int_t n ( numeric_cast<blas_int_t>( SO ? (*A).columns() : (*A).rows() ) );
395 blas_int_t lda ( numeric_cast<blas_int_t>( (*A).spacing() ) );
396 blas_int_t ldu ( numeric_cast<blas_int_t>( (*U).spacing() ) );
397 blas_int_t info( 0 );
398
399 const blas_int_t minimum( min( m, n ) );
400 const blas_int_t maximum( max( m, n ) );
401
402 blas_int_t lwork( ( jobz == 'O' )
403 ?( 2*minimum*minimum + 2*minimum + maximum + 2 )
404 :( 2*minimum + maximum + 2 ) );
405 const blas_int_t lrwork( ( jobz == 'O' )
406 ?( max( 5*minimum*minimum + 5*minimum,
407 2*maximum*minimum + 2*minimum*minimum + minimum ) )
408 :( 7*minimum ) );
409 const std::unique_ptr<CT[]> work( new CT[lwork] );
410 const std::unique_ptr<BT[]> rwork( new BT[lrwork] );
411 const std::unique_ptr<blas_int_t[]> iwork( new blas_int_t[8*minimum] );
412
413 gesdd( jobz, m, n, (*A).data(), lda, (*s).data(),
414 ( SO ? (*U).data() : nullptr ), ( SO ? ldu : 1 ),
415 ( SO ? nullptr : (*U).data() ), ( SO ? 1 : ldu ),
416 work.get(), lwork, rwork.get(), iwork.get(), &info );
417
418 BLAZE_INTERNAL_ASSERT( info >= 0, "Invalid argument for singular value decomposition" );
419
420 if( info > 0 ) {
421 BLAZE_THROW_LAPACK_ERROR( "Singular value decomposition failed" );
422 }
423}
425//*************************************************************************************************
426
427
428//*************************************************************************************************
517template< typename MT1 // Type of the matrix A
518 , bool SO // Storage order of all matrices
519 , typename MT2 // Type of the matrix U
520 , typename VT // Type of the vector s
521 , bool TF > // Transpose flag of the vector s
523 DenseVector<VT,TF>& s, char jobz )
524{
530
536
541
542 const size_t M( (*A).rows() );
543 const size_t N( (*A).columns() );
544 const size_t mindim( min( M, N ) );
545
546 if( jobz != 'O' && jobz != 'N' ) {
547 BLAZE_THROW_INVALID_ARGUMENT( "Invalid jobz argument provided" );
548 }
549
550 if( jobz == 'O' && M >= N ) {
551 BLAZE_THROW_INVALID_ARGUMENT( "Invalid input matrix provided" );
552 }
553
554 resize( *s, mindim, false );
555
556 if( jobz == 'O' ) {
557 resize( *U, M, M, false );
558 }
559
560 if( (*A).rows() == 0UL || (*A).columns() == 0UL ) {
561 return;
562 }
563
564 gesdd_backend( *A, *U, *s, jobz );
565}
566//*************************************************************************************************
567
568
569//*************************************************************************************************
588template< typename MT1 // Type of the matrix A
589 , bool SO // Storage order of all matrices
590 , typename VT // Type of the vector s
591 , bool TF // Transpose flag of the vector s
592 , typename MT2 > // Type of the matrix V
593inline auto gesdd_backend( DenseMatrix<MT1,SO>& A, DenseVector<VT,TF>& s,
594 DenseMatrix<MT2,SO>& V, char jobz )
595 -> DisableIf_t< IsComplex_v< ElementType_t<MT1> > >
596{
597 BLAZE_INTERNAL_ASSERT( jobz == 'O' || jobz == 'N', "Invalid jobz flag detected" );
598 BLAZE_INTERNAL_ASSERT( jobz == 'N' || isSquare( *V ), "Invalid matrix dimensions detected" );
599 BLAZE_INTERNAL_ASSERT( jobz == 'N' || (*V).rows() == (*A).columns(), "Invalid matrix dimensions detected" );
600 BLAZE_INTERNAL_ASSERT( (*s).size() == min( (*A).rows(), (*A).columns() ), "Invalid vector dimension detected" );
601
602 using ET = ElementType_t<MT1>;
603
604 blas_int_t m ( numeric_cast<blas_int_t>( SO ? (*A).rows() : (*A).columns() ) );
605 blas_int_t n ( numeric_cast<blas_int_t>( SO ? (*A).columns() : (*A).rows() ) );
606 blas_int_t lda ( numeric_cast<blas_int_t>( (*A).spacing() ) );
607 blas_int_t ldv ( numeric_cast<blas_int_t>( (*V).spacing() ) );
608 blas_int_t info( 0 );
609
610 const blas_int_t minimum( min( m, n ) );
611 const blas_int_t maximum( max( m, n ) );
612
613 blas_int_t lwork( ( jobz == 'O' )
614 ?( 3*minimum + max( maximum, 5*minimum*minimum + 4*maximum + 2 ) )
615 :( 3*minimum + max( maximum, 7*minimum ) + 2 ) );
616 const std::unique_ptr<ET[]> work( new ET[lwork] );
617 const std::unique_ptr<blas_int_t[]> iwork( new blas_int_t[8*minimum] );
618
619 gesdd( jobz, m, n, (*A).data(), lda, (*s).data(),
620 ( SO ? nullptr : (*V).data() ), ( SO ? 1 : ldv ),
621 ( SO ? (*V).data() : nullptr ), ( SO ? ldv : 1 ),
622 work.get(), lwork, iwork.get(), &info );
623
624 BLAZE_INTERNAL_ASSERT( info >= 0, "Invalid argument for singular value decomposition" );
625
626 if( info > 0 ) {
627 BLAZE_THROW_LAPACK_ERROR( "Singular value decomposition failed" );
628 }
629}
631//*************************************************************************************************
632
633
634//*************************************************************************************************
653template< typename MT1 // Type of the matrix A
654 , bool SO // Storage order of all matrices
655 , typename VT // Type of the vector s
656 , bool TF // Transpose flag of the vector s
657 , typename MT2 > // Type of the matrix V
658inline auto gesdd_backend( DenseMatrix<MT1,SO>& A, DenseVector<VT,TF>& s,
659 DenseMatrix<MT2,SO>& V, char jobz )
660 -> EnableIf_t< IsComplex_v< ElementType_t<MT1> > >
661{
662 BLAZE_INTERNAL_ASSERT( jobz == 'O' || jobz == 'N', "Invalid jobz flag detected" );
663 BLAZE_INTERNAL_ASSERT( jobz == 'N' || isSquare( *V ), "Invalid matrix dimensions detected" );
664 BLAZE_INTERNAL_ASSERT( jobz == 'N' || (*V).rows() == (*A).columns(), "Invalid matrix dimensions detected" );
665 BLAZE_INTERNAL_ASSERT( (*s).size() == min( (*A).rows(), (*A).columns() ), "Invalid vector dimension detected" );
666
667 using CT = ElementType_t<MT1>;
668 using BT = UnderlyingElement_t<CT>;
669
670 blas_int_t m ( numeric_cast<blas_int_t>( SO ? (*A).rows() : (*A).columns() ) );
671 blas_int_t n ( numeric_cast<blas_int_t>( SO ? (*A).columns() : (*A).rows() ) );
672 blas_int_t lda ( numeric_cast<blas_int_t>( (*A).spacing() ) );
673 blas_int_t ldv ( numeric_cast<blas_int_t>( (*V).spacing() ) );
674 blas_int_t info( 0 );
675
676 const blas_int_t minimum( min( m, n ) );
677 const blas_int_t maximum( max( m, n ) );
678
679 blas_int_t lwork( ( jobz == 'O' )
680 ?( 2*minimum*minimum + 2*minimum + maximum + 2 )
681 :( 2*minimum + maximum + 2 ) );
682 const blas_int_t lrwork( ( jobz == 'O' )
683 ?( max( 5*minimum*minimum + 5*minimum,
684 2*maximum*minimum + 2*minimum*minimum + minimum ) )
685 :( 7*minimum ) );
686 const std::unique_ptr<CT[]> work( new CT[lwork] );
687 const std::unique_ptr<BT[]> rwork( new BT[lrwork] );
688 const std::unique_ptr<blas_int_t[]> iwork( new blas_int_t[8*minimum] );
689
690 gesdd( jobz, m, n, (*A).data(), lda, (*s).data(),
691 ( SO ? nullptr : (*V).data() ), ( SO ? 1 : ldv ),
692 ( SO ? (*V).data() : nullptr ), ( SO ? ldv : 1 ),
693 work.get(), lwork, rwork.get(), iwork.get(), &info );
694
695 BLAZE_INTERNAL_ASSERT( info >= 0, "Invalid argument for singular value decomposition" );
696
697 if( info > 0 ) {
698 BLAZE_THROW_LAPACK_ERROR( "Singular value decomposition failed" );
699 }
700}
702//*************************************************************************************************
703
704
705//*************************************************************************************************
785template< typename MT1 // Type of the matrix A
786 , bool SO // Storage order of all matrices
787 , typename MT2 // Type of the matrix U
788 , typename VT // Type of the vector s
789 , bool TF > // Transpose flag of the vector s
791 DenseMatrix<MT2,SO>& V, char jobz )
792{
798
804
809
810 const size_t M( (*A).rows() );
811 const size_t N( (*A).columns() );
812 const size_t mindim( min( M, N ) );
813
814 if( jobz != 'O' && jobz != 'N' ) {
815 BLAZE_THROW_INVALID_ARGUMENT( "Invalid jobz argument provided" );
816 }
817
818 if( jobz == 'O' && M < N ) {
819 BLAZE_THROW_INVALID_ARGUMENT( "Invalid input matrix provided" );
820 }
821
822 resize( *s, mindim, false );
823
824 if( jobz == 'O' ) {
825 resize( *V, N, N, false );
826 }
827
828 if( (*A).rows() == 0UL || (*A).columns() == 0UL ) {
829 return;
830 }
831
832 gesdd_backend( *A, *s, *V, jobz );
833}
834//*************************************************************************************************
835
836
837//*************************************************************************************************
857template< typename MT1 // Type of the matrix A
858 , bool SO // Storage order of all matrices
859 , typename MT2 // Type of the matrix U
860 , typename VT // Type of the vector s
861 , bool TF // Transpose flag of the vector s
862 , typename MT3 > // Type of the matrix V
863inline auto gesdd_backend( DenseMatrix<MT1,SO>& A, DenseMatrix<MT2,SO>& U,
864 DenseVector<VT,TF>& s, DenseMatrix<MT3,SO>& V, char jobz )
865 -> DisableIf_t< IsComplex_v< ElementType_t<MT1> > >
866{
867 BLAZE_INTERNAL_ASSERT( jobz == 'A' || jobz == 'S' || jobz == 'N', "Invalid jobz flag detected" );
868 BLAZE_INTERNAL_ASSERT( jobz == 'N' || (*U).rows() == (*A).rows(), "Invalid matrix dimension detected" );
869 BLAZE_INTERNAL_ASSERT( jobz != 'A' || isSquare( *U ), "Invalid non-square matrix detected" );
870 BLAZE_INTERNAL_ASSERT( jobz != 'S' || (*U).columns() == min( (*A).rows(), (*A).columns() ), "Invalid matrix dimension detected" );
871 BLAZE_INTERNAL_ASSERT( jobz == 'N' || (*V).columns() == (*A).columns(), "Invalid matrix dimension detected" );
872 BLAZE_INTERNAL_ASSERT( jobz != 'A' || isSquare( *V ), "Invalid non-square matrix detected" );
873 BLAZE_INTERNAL_ASSERT( jobz != 'S' || (*V).rows() == min( (*A).rows(), (*A).columns() ), "Invalid matrix dimension detected" );
874 BLAZE_INTERNAL_ASSERT( (*s).size() == min( (*A).rows(), (*A).columns() ), "Invalid vector dimension detected" );
875
876 using ET = ElementType_t<MT1>;
877
878 blas_int_t m ( numeric_cast<blas_int_t>( SO ? (*A).rows() : (*A).columns() ) );
879 blas_int_t n ( numeric_cast<blas_int_t>( SO ? (*A).columns() : (*A).rows() ) );
880 blas_int_t lda ( numeric_cast<blas_int_t>( (*A).spacing() ) );
881 blas_int_t ldu ( numeric_cast<blas_int_t>( (*U).spacing() ) );
882 blas_int_t ldv ( numeric_cast<blas_int_t>( (*V).spacing() ) );
883 blas_int_t info( 0 );
884
885 const blas_int_t minimum( min( m, n ) );
886 const blas_int_t maximum( max( m, n ) );
887
888 blas_int_t lwork( 4*minimum*minimum + 6*minimum + maximum + 2 );
889 const std::unique_ptr<ET[]> work( new ET[lwork] );
890 const std::unique_ptr<blas_int_t[]> iwork( new blas_int_t[8*minimum] );
891
892 gesdd( jobz, m, n, (*A).data(), lda, (*s).data(),
893 ( SO ? (*U).data() : (*V).data() ), ( SO ? ldu : ldv ),
894 ( SO ? (*V).data() : (*U).data() ), ( SO ? ldv : ldu ),
895 work.get(), lwork, iwork.get(), &info );
896
897 BLAZE_INTERNAL_ASSERT( info >= 0, "Invalid argument for singular value decomposition" );
898
899 if( info > 0 ) {
900 BLAZE_THROW_LAPACK_ERROR( "Singular value decomposition failed" );
901 }
902}
904//*************************************************************************************************
905
906
907//*************************************************************************************************
927template< typename MT1 // Type of the matrix A
928 , bool SO // Storage order of all matrices
929 , typename MT2 // Type of the matrix U
930 , typename VT // Type of the vector s
931 , bool TF // Transpose flag of the vector s
932 , typename MT3 > // Type of the matrix V
933inline auto gesdd_backend( DenseMatrix<MT1,SO>& A, DenseMatrix<MT2,SO>& U,
934 DenseVector<VT,TF>& s, DenseMatrix<MT3,SO>& V, char jobz )
935 -> EnableIf_t< IsComplex_v< ElementType_t<MT1> > >
936{
937 BLAZE_INTERNAL_ASSERT( jobz == 'A' || jobz == 'S' || jobz == 'N', "Invalid jobz flag detected" );
938 BLAZE_INTERNAL_ASSERT( jobz == 'N' || (*U).rows() == (*A).rows(), "Invalid matrix dimension detected" );
939 BLAZE_INTERNAL_ASSERT( jobz != 'A' || isSquare( *U ), "Invalid non-square matrix detected" );
940 BLAZE_INTERNAL_ASSERT( jobz != 'S' || (*U).columns() == min( (*A).rows(), (*A).columns() ), "Invalid matrix dimension detected" );
941 BLAZE_INTERNAL_ASSERT( jobz == 'N' || (*V).columns() == (*A).columns(), "Invalid matrix dimension detected" );
942 BLAZE_INTERNAL_ASSERT( jobz != 'A' || isSquare( *V ), "Invalid non-square matrix detected" );
943 BLAZE_INTERNAL_ASSERT( jobz != 'S' || (*V).rows() == min( (*A).rows(), (*A).columns() ), "Invalid matrix dimension detected" );
944 BLAZE_INTERNAL_ASSERT( (*s).size() == min( (*A).rows(), (*A).columns() ), "Invalid vector dimension detected" );
945
946 using CT = ElementType_t<MT1>;
947 using BT = UnderlyingElement_t<CT>;
948
949 blas_int_t m ( numeric_cast<blas_int_t>( SO ? (*A).rows() : (*A).columns() ) );
950 blas_int_t n ( numeric_cast<blas_int_t>( SO ? (*A).columns() : (*A).rows() ) );
951 blas_int_t lda ( numeric_cast<blas_int_t>( (*A).spacing() ) );
952 blas_int_t ldu ( numeric_cast<blas_int_t>( (*U).spacing() ) );
953 blas_int_t ldv ( numeric_cast<blas_int_t>( (*V).spacing() ) );
954 blas_int_t info( 0 );
955
956 const blas_int_t minimum( min( m, n ) );
957 const blas_int_t maximum( max( m, n ) );
958
959 blas_int_t lwork( 4*minimum*minimum + 6*minimum + maximum + 2 );
960 const blas_int_t lrwork( max( 5*minimum*minimum + 5*minimum,
961 2*maximum*minimum + 2*minimum*minimum + minimum ) );
962 const std::unique_ptr<CT[]> work( new CT[lwork] );
963 const std::unique_ptr<BT[]> rwork( new BT[lrwork] );
964 const std::unique_ptr<blas_int_t[]> iwork( new blas_int_t[8*minimum] );
965
966 gesdd( jobz, m, n, (*A).data(), lda, (*s).data(),
967 ( SO ? (*U).data() : (*V).data() ), ( SO ? ldu : ldv ),
968 ( SO ? (*V).data() : (*U).data() ), ( SO ? ldv : ldu ),
969 work.get(), lwork, rwork.get(), iwork.get(), &info );
970
971 BLAZE_INTERNAL_ASSERT( info >= 0, "Invalid argument for singular value decomposition" );
972
973 if( info > 0 ) {
974 BLAZE_THROW_LAPACK_ERROR( "Singular value decomposition failed" );
975 }
976}
978//*************************************************************************************************
979
980
981//*************************************************************************************************
1083template< typename MT1 // Type of the matrix A
1084 , bool SO // Storage order of all matrices
1085 , typename MT2 // Type of the matrix U
1086 , typename VT // Type of the vector s
1087 , bool TF // Transpose flag of the vector s
1088 , typename MT3 > // Type of the matrix V
1090 DenseVector<VT,TF>& s, DenseMatrix<MT3,SO>& V, char jobz )
1091{
1097
1103
1108
1114
1115 const size_t M( (*A).rows() );
1116 const size_t N( (*A).columns() );
1117 const size_t mindim( min( M, N ) );
1118
1119 if( jobz != 'A' && jobz != 'S' && jobz != 'N' ) {
1120 BLAZE_THROW_INVALID_ARGUMENT( "Invalid jobz argument provided" );
1121 }
1122
1123 resize( *s, mindim, false );
1124
1125 if( jobz != 'N' ) {
1126 resize( *U, M, ( jobz == 'A' ? M : mindim ), false );
1127 resize( *V, ( jobz == 'A' ? N : mindim ), N, false );
1128 }
1129
1130 if( (*A).rows() == 0UL || (*A).columns() == 0UL ) {
1131 return;
1132 }
1133
1134 gesdd_backend( *A, *U, *s, *V, jobz );
1135}
1136//*************************************************************************************************
1137
1138} // namespace blaze
1139
1140#endif
Constraint on the data type.
Header file for auxiliary alias declarations.
typename T::ElementType ElementType_t
Alias declaration for nested ElementType type definitions.
Definition: Aliases.h:190
Header file for run time assertion macros.
Constraint on the data type.
Constraint on the data type.
Header file for the EnableIf class template.
Header file for the IsComplex type trait.
Constraint on the data type.
Cast operators for numeric types.
Header file for the UnderlyingElement type trait.
Header file for the CLAPACK gesdd wrapper functions.
Base class for dense matrices.
Definition: DenseMatrix.h:82
Base class for N-dimensional dense vectors.
Definition: DenseVector.h:77
Constraint on the data type.
Header file for the DenseMatrix base class.
Header file for the DenseVector base class.
decltype(auto) min(const DenseMatrix< MT1, SO1 > &lhs, const DenseMatrix< MT2, SO2 > &rhs)
Computes the componentwise minimum of the dense matrices lhs and rhs.
Definition: DMatDMatMapExpr.h:1339
decltype(auto) max(const DenseMatrix< MT1, SO1 > &lhs, const DenseMatrix< MT2, SO2 > &rhs)
Computes the componentwise maximum of the dense matrices lhs and rhs.
Definition: DMatDMatMapExpr.h:1375
void gesdd(DenseMatrix< MT1, SO > &A, DenseMatrix< MT2, SO > &U, DenseVector< VT, TF > &s, DenseMatrix< MT3, SO > &V, char jobz)
LAPACK kernel for the singular value decomposition (SVD) of the given dense general matrix.
Definition: gesdd.h:1089
#define BLAZE_CONSTRAINT_MUST_BE_BLAS_COMPATIBLE_TYPE(T)
Constraint on the data type.
Definition: BLASCompatible.h:61
#define BLAZE_CONSTRAINT_MUST_BE_CONTIGUOUS_TYPE(T)
Constraint on the data type.
Definition: Contiguous.h:61
#define BLAZE_CONSTRAINT_MUST_NOT_BE_COMPUTATION_TYPE(T)
Constraint on the data type.
Definition: Computation.h:81
#define BLAZE_CONSTRAINT_MUST_NOT_BE_ADAPTOR_TYPE(T)
Constraint on the data type.
Definition: Adaptor.h:81
#define BLAZE_CONSTRAINT_MUST_HAVE_MUTABLE_DATA_ACCESS(T)
Constraint on the data type.
Definition: MutableDataAccess.h:61
int32_t blas_int_t
Signed integer type used in the BLAS/LAPACK wrapper functions.
Definition: Types.h:64
#define BLAZE_THROW_LAPACK_ERROR(MESSAGE)
Macro for the emission of an exception on detection of a LAPACK error.
Definition: Exception.h:146
void resize(Matrix< MT, SO > &matrix, size_t rows, size_t columns, bool preserve=true)
Changing the size of the matrix.
Definition: Matrix.h:1108
bool isSquare(const Matrix< MT, SO > &matrix) noexcept
Checks if the given matrix is a square matrix.
Definition: Matrix.h:1383
#define BLAZE_INTERNAL_ASSERT(expr, msg)
Run time assertion macro for internal checks.
Definition: Assert.h:101
#define BLAZE_THROW_INVALID_ARGUMENT(MESSAGE)
Macro for the emission of a std::invalid_argument exception.
Definition: Exception.h:235
Header file for the exception macros of the math module.
Header file for the generic max algorithm.
Header file for the generic min algorithm.