QDWHInfo QDWH( Matrix<F>& A, const QDWHCtrl& ctrl ) { EL_DEBUG_CSE typedef Base<F> Real; const Real twoEst = TwoNormEstimate( A ); A *= 1/twoEst; // The one-norm of the inverse can be replaced with an estimate which is // a few times cheaper, e.g., via Higham and Tisseur's block algorithm // from "A Block Algorithm for Matrix 1-Norm Estimation, with an Application // to 1-Norm Pseudospectra". Real sMinUpper; Matrix<F> Y( A ); if( A.Height() > A.Width() ) { qr::ExplicitTriang( Y ); try { TriangularInverse( UPPER, NON_UNIT, Y ); sMinUpper = Real(1) / OneNorm( Y ); } catch( SingularMatrixException& e ) { sMinUpper = 0; } } else { try { Inverse( Y ); sMinUpper = Real(1) / OneNorm( Y ); } catch( SingularMatrixException& e ) { sMinUpper = 0; } } return QDWHInner( A, sMinUpper, ctrl ); }
inline void CholeskyUVar2( Matrix<F>& A ) { #ifndef RELEASE PushCallStack("hpd_inverse::CholeskyUVar2"); if( A.Height() != A.Width() ) throw std::logic_error("Nonsquare matrices cannot be triangular"); #endif // Matrix views Matrix<F> ATL, ATR, A00, A01, A02, ABL, ABR, A10, A11, A12, A20, A21, A22; // Start the algorithm PartitionDownDiagonal ( A, ATL, ATR, ABL, ABR, 0 ); while( ATL.Height() < A.Height() ) { RepartitionDownDiagonal ( ATL, /**/ ATR, A00, /**/ A01, A02, /*************/ /******************/ /**/ A10, /**/ A11, A12, ABL, /**/ ABR, A20, /**/ A21, A22 ); //--------------------------------------------------------------------// Cholesky( UPPER, A11 ); Trsm( RIGHT, UPPER, NORMAL, NON_UNIT, F(1), A11, A01 ); Trsm( LEFT, UPPER, ADJOINT, NON_UNIT, F(1), A11, A12 ); Herk( UPPER, NORMAL, F(1), A01, F(1), A00 ); Gemm( NORMAL, NORMAL, F(-1), A01, A12, F(1), A02 ); Herk( UPPER, ADJOINT, F(-1), A12, F(1), A22 ); Trsm( RIGHT, UPPER, ADJOINT, NON_UNIT, F(1), A11, A01 ); Trsm( LEFT, UPPER, NORMAL, NON_UNIT, F(-1), A11, A12 ); TriangularInverse( UPPER, NON_UNIT, A11 ); Trtrmm( ADJOINT, UPPER, A11 ); //--------------------------------------------------------------------// SlidePartitionDownDiagonal ( ATL, /**/ ATR, A00, A01, /**/ A02, /**/ A10, A11, /**/ A12, /*************/ /******************/ ABL, /**/ ABR, A20, A21, /**/ A22 ); } #ifndef RELEASE PopCallStack(); #endif }
Int ADMM ( const Matrix<Real>& A, const Matrix<Real>& b, const Matrix<Real>& c, Matrix<Real>& z, const ADMMCtrl<Real>& ctrl ) { EL_DEBUG_CSE // Cache a custom partially-pivoted LU factorization of // | rho*I A^H | = | B11 B12 | // | A 0 | | B21 B22 | // by (justifiably) avoiding pivoting in the first n steps of // the factorization, so that // [I,rho*I] = lu(rho*I). // The factorization would then proceed with // B21 := B21 U11^{-1} = A (rho*I)^{-1} = A/rho // B12 := L11^{-1} B12 = I A^H = A^H. // The Schur complement would then be // B22 := B22 - B21 B12 = 0 - (A*A^H)/rho. // We then factor said matrix with LU with partial pivoting and // swap the necessary rows of B21 in order to implicitly commute // the row pivots with the Gauss transforms in the manner standard // for GEPP. Unless A A' is singular, pivoting should not be needed, // as Cholesky factorization of the negative matrix should be valid. // // The result is the factorization // | I 0 | | rho*I A^H | = | I 0 | | rho*I U12 |, // | 0 P22 | | A 0 | | L21 L22 | | 0 U22 | // where [L22,U22] are stored within B22. Matrix<Real> U12, L21, B22, bPiv; Adjoint( A, U12 ); L21 = A; L21 *= 1/ctrl.rho; Herk( LOWER, NORMAL, -1/ctrl.rho, A, B22 ); MakeHermitian( LOWER, B22 ); // TODO: Replace with sparse-direct Cholesky version? Permutation P2; LU( B22, P2 ); P2.PermuteRows( L21 ); bPiv = b; P2.PermuteRows( bPiv ); // Possibly form the inverse of L22 U22 Matrix<Real> X22; if( ctrl.inv ) { X22 = B22; MakeTrapezoidal( LOWER, X22 ); FillDiagonal( X22, Real(1) ); TriangularInverse( LOWER, UNIT, X22 ); Trsm( LEFT, UPPER, NORMAL, NON_UNIT, Real(1), B22, X22 ); } Int numIter=0; const Int m = A.Height(); const Int n = A.Width(); Matrix<Real> g, xTmp, y, t; Zeros( g, m+n, 1 ); PartitionDown( g, xTmp, y, n ); Matrix<Real> x, u, zOld, xHat; Zeros( z, n, 1 ); Zeros( u, n, 1 ); Zeros( t, n, 1 ); while( numIter < ctrl.maxIter ) { zOld = z; // Find x from // | rho*I A^H | | x | = | rho*(z-u)-c | // | A 0 | | y | | b | // via our cached custom factorization: // // |x| = inv(U) inv(L) P' |rho*(z-u)-c| // |y| |b | // = |rho*I U12|^{-1} |I 0 | |I 0 | |rho*(z-u)-c| // = |0 U22| |L21 L22| |0 P22'| |b | // = " " |rho*(z-u)-c| // | P22' b | xTmp = z; xTmp -= u; xTmp *= ctrl.rho; xTmp -= c; y = bPiv; Gemv( NORMAL, Real(-1), L21, xTmp, Real(1), y ); if( ctrl.inv ) { Gemv( NORMAL, Real(1), X22, y, t ); y = t; } else { Trsv( LOWER, NORMAL, UNIT, B22, y ); Trsv( UPPER, NORMAL, NON_UNIT, B22, y ); } Gemv( NORMAL, Real(-1), U12, y, Real(1), xTmp ); xTmp *= 1/ctrl.rho; // xHat := alpha*x + (1-alpha)*zOld xHat = xTmp; xHat *= ctrl.alpha; Axpy( 1-ctrl.alpha, zOld, xHat ); // z := pos(xHat+u) z = xHat; z += u; LowerClip( z, Real(0) ); // u := u + (xHat-z) u += xHat; u -= z; const Real objective = Dot( c, xTmp ); // rNorm := || x - z ||_2 t = xTmp; t -= z; const Real rNorm = FrobeniusNorm( t ); // sNorm := |rho| || z - zOld ||_2 t = z; t -= zOld; const Real sNorm = Abs(ctrl.rho)*FrobeniusNorm( t ); const Real epsPri = Sqrt(Real(n))*ctrl.absTol + ctrl.relTol*Max(FrobeniusNorm(xTmp),FrobeniusNorm(z)); const Real epsDual = Sqrt(Real(n))*ctrl.absTol + ctrl.relTol*Abs(ctrl.rho)*FrobeniusNorm(u); if( ctrl.print ) { t = xTmp; LowerClip( t, Real(0) ); t -= xTmp; const Real clipDist = FrobeniusNorm( t ); cout << numIter << ": " << "||x-z||_2=" << rNorm << ", " << "epsPri=" << epsPri << ", " << "|rho| ||z-zOld||_2=" << sNorm << ", " << "epsDual=" << epsDual << ", " << "||x-Pos(x)||_2=" << clipDist << ", " << "c'x=" << objective << endl; } if( rNorm < epsPri && sNorm < epsDual ) break; ++numIter; } if( ctrl.maxIter == numIter ) cout << "ADMM failed to converge" << endl; x = xTmp; return numIter; }
inline void Inverse( DistMatrix<F>& A ) { #ifndef RELEASE PushCallStack("Inverse"); if( A.Height() != A.Width() ) throw std::logic_error("Cannot invert non-square matrices"); #endif const Grid& g = A.Grid(); DistMatrix<int,VC,STAR> p( g ); LU( A, p ); TriangularInverse( UPPER, NON_UNIT, A ); // Solve inv(A) L = inv(U) for inv(A) DistMatrix<F> ATL(g), ATR(g), ABL(g), ABR(g); DistMatrix<F> A00(g), A01(g), A02(g), A10(g), A11(g), A12(g), A20(g), A21(g), A22(g); DistMatrix<F> A1(g), A2(g); DistMatrix<F,VC, STAR> A1_VC_STAR(g); DistMatrix<F,STAR,STAR> L11_STAR_STAR(g); DistMatrix<F,VR, STAR> L21_VR_STAR(g); DistMatrix<F,STAR,MR > L21Trans_STAR_MR(g); DistMatrix<F,MC, STAR> Z1(g); PartitionUpDiagonal ( A, ATL, ATR, ABL, ABR, 0 ); while( ABR.Height() < A.Height() ) { RepartitionUpDiagonal ( ATL, /**/ ATR, A00, A01, /**/ A02, /**/ A10, A11, /**/ A12, /*************/ /******************/ ABL, /**/ ABR, A20, A21, /**/ A22 ); View( A1, A, 0, A00.Width(), A.Height(), A01.Width() ); View( A2, A, 0, A00.Width()+A01.Width(), A.Height(), A02.Width() ); L21_VR_STAR.AlignWith( A2 ); L21Trans_STAR_MR.AlignWith( A2 ); Z1.AlignWith( A01 ); Z1.ResizeTo( A.Height(), A01.Width() ); //--------------------------------------------------------------------// // Copy out L1 L11_STAR_STAR = A11; L21_VR_STAR = A21; L21Trans_STAR_MR.TransposeFrom( L21_VR_STAR ); // Zero the strictly lower triangular portion of A1 MakeTrapezoidal( LEFT, UPPER, 0, A11 ); Zero( A21 ); // Perform the lazy update of A1 internal::LocalGemm ( NORMAL, TRANSPOSE, F(-1), A2, L21Trans_STAR_MR, F(0), Z1 ); A1.SumScatterUpdate( F(1), Z1 ); // Solve against this diagonal block of L11 A1_VC_STAR = A1; internal::LocalTrsm ( RIGHT, LOWER, NORMAL, UNIT, F(1), L11_STAR_STAR, A1_VC_STAR ); A1 = A1_VC_STAR; //--------------------------------------------------------------------// Z1.FreeAlignments(); L21Trans_STAR_MR.FreeAlignments(); L21_VR_STAR.FreeAlignments(); SlidePartitionUpDiagonal ( ATL, /**/ ATR, A00, /**/ A01, A02, /*************/ /*******************/ /**/ A10, /**/ A11, A12, ABL, /**/ ABR, A20, /**/ A21, A22 ); } // inv(A) := inv(A) P ApplyInverseColumnPivots( A, p ); #ifndef RELEASE PopCallStack(); #endif }
inline void Inverse( Matrix<F>& A ) { #ifndef RELEASE PushCallStack("Inverse"); if( A.Height() != A.Width() ) throw std::logic_error("Cannot invert non-square matrices"); #endif Matrix<int> p; LU( A, p ); TriangularInverse( UPPER, NON_UNIT, A ); // Solve inv(A) L = inv(U) for inv(A) Matrix<F> ATL, ATR, ABL, ABR; Matrix<F> A00, A01, A02, A10, A11, A12, A20, A21, A22; Matrix<F> A1, A2; Matrix<F> L11, L21; PartitionUpDiagonal ( A, ATL, ATR, ABL, ABR, 0 ); while( ABR.Height() < A.Height() ) { RepartitionUpDiagonal ( ATL, /**/ ATR, A00, A01, /**/ A02, /**/ A10, A11, /**/ A12, /*************/ /******************/ ABL, /**/ ABR, A20, A21, /**/ A22 ); View( A1, A, 0, A00.Width(), A.Height(), A01.Width() ); View( A2, A, 0, A00.Width()+A01.Width(), A.Height(), A02.Width() ); //--------------------------------------------------------------------// // Copy out L1 L11 = A11; L21 = A21; // Zero the strictly lower triangular portion of A1 MakeTrapezoidal( LEFT, UPPER, 0, A11 ); Zero( A21 ); // Perform the lazy update of A1 Gemm( NORMAL, NORMAL, F(-1), A2, L21, F(1), A1 ); // Solve against this diagonal block of L11 Trsm( RIGHT, LOWER, NORMAL, UNIT, F(1), L11, A1 ); //--------------------------------------------------------------------// SlidePartitionUpDiagonal ( ATL, /**/ ATR, A00, /**/ A01, A02, /*************/ /*******************/ /**/ A10, /**/ A11, A12, ABL, /**/ ABR, A20, /**/ A21, A22 ); } // inv(A) := inv(A) P ApplyInverseColumnPivots( A, p ); #ifndef RELEASE PopCallStack(); #endif }