void ExplicitUnitary ( ElementalMatrix<F>& APre, bool thinQR, const QRCtrl<Base<F>>& ctrl ) { DEBUG_CSE DistMatrixReadWriteProxy<F,F,MC,MR> AProx( APre ); auto& A = AProx.Get(); const Grid& g = A.Grid(); DistMatrix<F,MD,STAR> householderScalars(g); DistMatrix<Base<F>,MD,STAR> signature(g); if( ctrl.colPiv ) { DistPermutation Omega(g); QR( A, householderScalars, signature, Omega, ctrl ); } else QR( A, householderScalars, signature ); if( thinQR ) { A.Resize( A.Height(), householderScalars.Height() ); ExpandPackedReflectors ( LOWER, VERTICAL, CONJUGATED, 0, A, householderScalars ); DiagonalScale( RIGHT, NORMAL, signature, A ); } else { auto ACopy = A; // TODO: Use an extension of ExpandPackedReflectors to make this faster Identity( A, A.Height(), A.Height() ); qr::ApplyQ( LEFT, NORMAL, ACopy, householderScalars, signature, A ); } }
void Explicit( AbstractDistMatrix<F>& L, AbstractDistMatrix<F>& APre ) { EL_DEBUG_CSE const Grid& g = APre.Grid(); DistMatrixReadWriteProxy<F,F,MC,MR> AProx( APre ); auto& A = AProx.Get(); DistMatrix<F,MD,STAR> householderScalars(g); DistMatrix<Base<F>,MD,STAR> signature(g); LQ( A, householderScalars, signature ); const Int m = A.Height(); const Int n = A.Width(); const Int minDim = Min(m,n); auto AL = A( IR(0,m), IR(0,minDim) ); Copy( AL, L ); MakeTrapezoidal( LOWER, L ); // TODO: Replace this with an in-place expansion of Q DistMatrix<F> Q(g); Identity( Q, A.Height(), A.Width() ); lq::ApplyQ( RIGHT, NORMAL, A, householderScalars, signature, Q ); Copy( Q, APre ); }
void ExplicitTriang( AbstractDistMatrix<F>& A ) { EL_DEBUG_CSE const Grid& g = A.Grid(); DistMatrix<F,MD,STAR> householderScalars(g); DistMatrix<Base<F>,MD,STAR> signature(g); LQ( A, householderScalars, signature ); const Int m = A.Height(); const Int n = A.Width(); const Int minDim = Min(m,n); A.Resize( m, minDim ); MakeTrapezoidal( LOWER, A ); }
void ExplicitTriang( ElementalMatrix<F>& A, const QRCtrl<Base<F>>& ctrl ) { DEBUG_CSE DistMatrix<F,MD,STAR> householderScalars(A.Grid()); DistMatrix<Base<F>,MD,STAR> signature(A.Grid()); if( ctrl.colPiv ) { DistPermutation Omega(A.Grid()); BusingerGolub( A, householderScalars, signature, Omega, ctrl ); } else Householder( A, householderScalars, signature ); A.Resize( householderScalars.Height(), A.Width() ); MakeTrapezoidal( UPPER, A ); }
void Householder ( AbstractDistMatrix<F>& APre, AbstractDistMatrix<F>& householderScalarsPre, AbstractDistMatrix<Base<F>>& signaturePre ) { DEBUG_CSE DEBUG_ONLY(AssertSameGrids( APre, householderScalarsPre, signaturePre )) const Int m = APre.Height(); const Int n = APre.Width(); const Int minDim = Min(m,n); const Int iOff = m-minDim; const Int jOff = n-minDim; DistMatrixReadWriteProxy<F,F,MC,MR> AProx( APre ); DistMatrixWriteProxy<F,F,MD,STAR> householderScalarsProx( householderScalarsPre ); DistMatrixWriteProxy<Base<F>,Base<F>,MD,STAR> signatureProx( signaturePre ); auto& A = AProx.Get(); auto& householderScalars = householderScalarsProx.Get(); auto& signature = signatureProx.Get(); householderScalars.Resize( minDim, 1 ); signature.Resize( minDim, 1 ); const Int bsize = Blocksize(); const Int kLast = LastOffset( minDim, bsize ); for( Int k=kLast; k>=0; k-=bsize ) { const Int nb = Min(bsize,minDim-k); const Int ki = k + iOff; const Int kj = k + jOff; const Range<Int> ind0Vert( 0, ki ), ind1( k, k+nb ), ind1Vert( ki, ki+nb ), indL( 0, kj+nb ); auto A0L = A( ind0Vert, indL ); auto A1L = A( ind1Vert, indL ); auto householderScalars1 = householderScalars( ind1, ALL ); auto sig1 = signature( ind1, ALL ); PanelHouseholder( A1L, householderScalars1, sig1 ); ApplyQ( RIGHT, ADJOINT, A1L, householderScalars1, sig1, A0L ); } }
void ExplicitUnitary( AbstractDistMatrix<F>& APre ) { EL_DEBUG_CSE const Grid& g = APre.Grid(); DistMatrixReadWriteProxy<F,F,MC,MR> AProx( APre ); auto& A = AProx.Get(); DistMatrix<F,MD,STAR> householderScalars(g); DistMatrix<Base<F>,MD,STAR> signature(g); LQ( A, householderScalars, signature ); // TODO: Replace this with an in-place expansion of Q DistMatrix<F> Q(g); Q.AlignWith( A ); Identity( Q, A.Height(), A.Width() ); lq::ApplyQ( RIGHT, NORMAL, A, householderScalars, signature, Q ); Copy( Q, APre ); }
void Explicit ( ElementalMatrix<F>& APre, ElementalMatrix<F>& R, ElementalMatrix<Int>& OmegaFull, bool thinQR, const QRCtrl<Base<F>>& ctrl ) { DEBUG_CSE DistMatrixReadWriteProxy<F,F,MC,MR> AProx( APre ); auto& A = AProx.Get(); const Grid& g = A.Grid(); DistMatrix<F,MD,STAR> householderScalars(g); DistMatrix<Base<F>,MD,STAR> signature(g); DistPermutation Omega(g); QR( A, householderScalars, signature, Omega, ctrl ); const Int m = A.Height(); const Int n = A.Width(); const Int numIts = householderScalars.Height(); auto AT = A( IR(0,numIts), IR(0,n) ); Copy( AT, R ); MakeTrapezoidal( UPPER, R ); if( thinQR ) { A.Resize( m, numIts ); ExpandPackedReflectors ( LOWER, VERTICAL, CONJUGATED, 0, A, householderScalars ); DiagonalScale( RIGHT, NORMAL, signature, A ); } else { auto ACopy = A; // TODO: Use an extension of ExpandPackedReflectors to make this faster Identity( A, A.Height(), A.Height() ); qr::ApplyQ( LEFT, NORMAL, ACopy, householderScalars, signature, A ); } Omega.ExplicitMatrix( OmegaFull ); }
void Skeleton ( const AbstractDistMatrix<F>& APre, DistPermutation& PR, DistPermutation& PC, AbstractDistMatrix<F>& Z, const QRCtrl<Base<F>>& ctrl ) { DEBUG_CSE DistMatrixReadProxy<F,F,MC,MR> AProx( APre ); auto& A = AProx.GetLocked(); const Grid& g = A.Grid(); DistMatrix<F> AAdj(g); Adjoint( A, AAdj ); // Find the row permutation DistMatrix<F> B(AAdj); DistMatrix<F,MD,STAR> householderScalars(g); DistMatrix<Base<F>,MD,STAR> signature(g); QR( B, householderScalars, signature, PR, ctrl ); const Int numSteps = householderScalars.Height(); B.Resize( B.Height(), numSteps ); // Form K' := (A pinv(AR))' = pinv(AR') A' DistMatrix<F> KAdj(g); qr::SolveAfter( NORMAL, B, householderScalars, signature, AAdj, KAdj ); // Form K := (K')' DistMatrix<F> K(g); Adjoint( KAdj, K ); // Find the column permutation (force the same number of steps) B = A; auto secondCtrl = ctrl; secondCtrl.adaptive = false; secondCtrl.boundRank = true; secondCtrl.maxRank = numSteps; QR( B, householderScalars, signature, PC, secondCtrl ); // Form Z := pinv(AC) K = pinv(AC) (A pinv(AR)) B.Resize( B.Height(), numSteps ); qr::SolveAfter( NORMAL, B, householderScalars, signature, K, Z ); }
inline void BusingerGolub ( AbstractDistMatrix<F>& APre, DistPermutation& Omega, AbstractDistMatrix<F>& Z, const QRCtrl<Base<F>>& ctrl ) { EL_DEBUG_CSE typedef Base<F> Real; DistMatrixReadWriteProxy<F,F,MC,MR> AProx( APre ); auto& A = AProx.Get(); auto ctrlCopy = ctrl; const Int m = A.Height(); const Int n = A.Width(); const Real eps = limits::Epsilon<Real>(); // Demand that we will be able to apply inv(R_L) to R_R by ensuring that // the minimum singular value is sufficiently (relatively) large ctrlCopy.adaptive = true; if( ctrl.boundRank ) { ctrlCopy.tol = Max(ctrl.tol,eps*ctrl.maxRank); } else { ctrlCopy.tol = Max(ctrl.tol,eps*Min(m,n)); } // Perform an adaptive pivoted QR factorization DistMatrix<F,MD,STAR> householderScalars(A.Grid()); DistMatrix<Base<F>,MD,STAR> signature(A.Grid()); QR( A, householderScalars, signature, Omega, ctrlCopy ); const Int numSteps = householderScalars.Height(); auto RL = A( IR(0,numSteps), IR(0,numSteps) ); auto RR = A( IR(0,numSteps), IR(numSteps,n) ); Copy( RR, Z ); Trsm( LEFT, UPPER, NORMAL, NON_UNIT, F(1), RL, Z ); }
void LowerBlocked ( AbstractDistMatrix<F>& APre, AbstractDistMatrix<F>& householderScalarsPre ) { DEBUG_CSE DEBUG_ONLY(AssertSameGrids( APre, householderScalarsPre )) DistMatrixReadWriteProxy<F,F,MC,MR> AProx( APre ); DistMatrixWriteProxy<F,F,STAR,STAR> householderScalarsProx( householderScalarsPre ); auto& A = AProx.Get(); auto& householderScalars = householderScalarsProx.Get(); const Grid& g = A.Grid(); const Int n = A.Height(); householderScalars.Resize( Max(n-1,0), 1 ); DistMatrix<F,MC,STAR> UB1_MC_STAR(g), V21_MC_STAR(g); DistMatrix<F,MR,STAR> V01_MR_STAR(g), VB1_MR_STAR(g), UB1_MR_STAR(g); DistMatrix<F,STAR,STAR> G11_STAR_STAR(g); const Int bsize = Blocksize(); for( Int k=0; k<n-1; k+=bsize ) { const Int nb = Min(bsize,n-1-k); const Range<Int> ind0( 0, k ), ind1( k, k+nb ), indB( k, n ), indR( k, n ), ind2( k+nb, n ); auto ABR = A( indB, indR ); auto A22 = A( ind2, ind2 ); auto householderScalars1 = householderScalars( ind1, ALL ); UB1_MC_STAR.AlignWith( ABR ); UB1_MR_STAR.AlignWith( ABR ); VB1_MR_STAR.AlignWith( ABR ); UB1_MC_STAR.Resize( n-k, nb ); UB1_MR_STAR.Resize( n-k, nb ); VB1_MR_STAR.Resize( n-k, nb ); G11_STAR_STAR.Resize( nb, nb ); hessenberg::LowerPanel ( ABR, householderScalars1, UB1_MC_STAR, UB1_MR_STAR, VB1_MR_STAR, G11_STAR_STAR ); auto AB0 = A( indB, ind0 ); auto A2R = A( ind2, indR ); auto U21_MC_STAR = UB1_MC_STAR( IR(nb,END), ALL ); // AB0 := AB0 - (UB1 inv(G11)^H UB1^H AB0) // = AB0 - (UB1 ((AB0^H UB1) inv(G11))^H) // ------------------------------------------- V01_MR_STAR.AlignWith( AB0 ); Zeros( V01_MR_STAR, k, nb ); LocalGemm( ADJOINT, NORMAL, F(1), AB0, UB1_MC_STAR, F(0), V01_MR_STAR ); El::AllReduce( V01_MR_STAR, AB0.ColComm() ); LocalTrsm ( RIGHT, UPPER, NORMAL, NON_UNIT, F(1), G11_STAR_STAR, V01_MR_STAR ); LocalGemm ( NORMAL, ADJOINT, F(-1), UB1_MC_STAR, V01_MR_STAR, F(1), AB0 ); // A2R := (A2R - U21 inv(G11)^H VB1^H)(I - UB1 inv(G11) UB1^H) // ----------------------------------------------------------- // A2R := A2R - U21 inv(G11)^H VB1^H // (note: VB1 is overwritten) LocalTrsm ( RIGHT, UPPER, NORMAL, NON_UNIT, F(1), G11_STAR_STAR, VB1_MR_STAR ); LocalGemm ( NORMAL, ADJOINT, F(-1), U21_MC_STAR, VB1_MR_STAR, F(1), A2R ); // A2R := A2R - ((A2R UB1) inv(G11)) UB1^H V21_MC_STAR.AlignWith( A2R ); Zeros( V21_MC_STAR, A2R.Height(), nb ); LocalGemm( NORMAL, NORMAL, F(1), A2R, UB1_MR_STAR, F(0), V21_MC_STAR ); El::AllReduce( V21_MC_STAR, A2R.RowComm() ); LocalTrsm ( RIGHT, UPPER, NORMAL, NON_UNIT, F(1), G11_STAR_STAR, V21_MC_STAR ); LocalGemm ( NORMAL, ADJOINT, F(-1), V21_MC_STAR, UB1_MR_STAR, F(1), A2R ); } }