int Trsm_unb_var2( FLA_Obj L, FLA_Obj B ) { FLA_Obj LTL, LTR, L00, l01, L02, LBL, LBR, l10t, lambda11, l12t, L20, l21, L22; FLA_Obj BT, B0, BB, b1t, B2; FLA_Part_2x2( L, <L, <R, &LBL, &LBR, 0, 0, FLA_TL ); FLA_Part_2x1( B, &BT, &BB, 0, FLA_TOP ); while ( FLA_Obj_length( LTL ) < FLA_Obj_length( L ) ){ FLA_Repart_2x2_to_3x3( LTL, /**/ LTR, &L00, /**/ &l01, &L02, /* ************* */ /* *************************** */ &l10t, /**/ &lambda11, &l12t, LBL, /**/ LBR, &L20, /**/ &l21, &L22, 1, 1, FLA_BR ); FLA_Repart_2x1_to_3x1( BT, &B0, /* ** */ /* *** */ &b1t, BB, &B2, 1, FLA_BOTTOM ); /*------------------------------------------------------------*/ /* b1t = b1t / lambda11 */ FLA_Inv_scal( lambda11, b1t ); /* B2 = B2 - l21 * b1t */ FLA_Ger( FLA_MINUS_ONE, l21, b1t, B2 ); /*------------------------------------------------------------*/ FLA_Cont_with_3x3_to_2x2( <L, /**/ <R, L00, l01, /**/ L02, l10t, lambda11, /**/ l12t, /* ************** */ /* ************************* */ &LBL, /**/ &LBR, L20, l21, /**/ L22, FLA_TL ); FLA_Cont_with_3x1_to_2x1( &BT, B0, b1t, /* ** */ /* *** */ &BB, B2, FLA_TOP ); } return FLA_SUCCESS; }
FLA_Error LU_unb_var4( FLA_Obj A ) { FLA_Obj ATL, ATR, A00, a01, A02, ABL, ABR, a10t, alpha11, a12t, A20, a21, A22; FLA_Part_2x2( A, &ATL, &ATR, &ABL, &ABR, 0, 0, FLA_TL ); while ( FLA_Obj_length( ATL ) < FLA_Obj_length( A ) ){ FLA_Repart_2x2_to_3x3( ATL, /**/ ATR, &A00, /**/ &a01, &A02, /* ************* */ /* *************************** */ &a10t, /**/ &alpha11, &a12t, ABL, /**/ ABR, &A20, /**/ &a21, &A22, 1, 1, FLA_BR ); /*------------------------------------------------------------*/ /* alpha11 = alpha11 - a10t * a01; */ FLA_Dots( FLA_MINUS_ONE, a10t, a01, FLA_ONE, alpha11 ); /* a12t = a12t - a10t * A02; */ FLA_Gemv( FLA_TRANSPOSE, FLA_MINUS_ONE, A02, a10t, FLA_ONE, a12t ); /* a21 = a21 - A20 * a01; */ FLA_Gemv( FLA_NO_TRANSPOSE, FLA_MINUS_ONE, A20, a01, FLA_ONE, a21 ); /* a21 = a21 / alpha11; */ FLA_Inv_scal( alpha11, a21 ); /*------------------------------------------------------------*/ FLA_Cont_with_3x3_to_2x2( &ATL, /**/ &ATR, A00, a01, /**/ A02, a10t, alpha11, /**/ a12t, /* ************** */ /* ************************* */ &ABL, /**/ &ABR, A20, a21, /**/ A22, FLA_TL ); } return FLA_SUCCESS; }
FLA_Error FLA_Svd_uv_unb_var1( dim_t n_iter_max, FLA_Obj A, FLA_Obj s, FLA_Obj U, FLA_Obj V, dim_t k_accum, dim_t b_alg ) { FLA_Error r_val = FLA_SUCCESS; FLA_Datatype dt; FLA_Datatype dt_real; FLA_Datatype dt_comp; FLA_Obj scale, T, S, rL, rR, d, e, G, H; dim_t m_A, n_A; dim_t min_m_n; dim_t n_GH; double crossover_ratio = 17.0 / 9.0; n_GH = k_accum; m_A = FLA_Obj_length( A ); n_A = FLA_Obj_width( A ); min_m_n = FLA_Obj_min_dim( A ); dt = FLA_Obj_datatype( A ); dt_real = FLA_Obj_datatype_proj_to_real( A ); dt_comp = FLA_Obj_datatype_proj_to_complex( A ); // Create matrices to hold block Householder transformations. FLA_Bidiag_UT_create_T( A, &T, &S ); // Create vectors to hold the realifying scalars. FLA_Obj_create( dt, min_m_n, 1, 0, 0, &rL ); FLA_Obj_create( dt, min_m_n, 1, 0, 0, &rR ); // Create vectors to hold the diagonal and sub-diagonal. FLA_Obj_create( dt_real, min_m_n, 1, 0, 0, &d ); FLA_Obj_create( dt_real, min_m_n-1, 1, 0, 0, &e ); // Create matrices to hold the left and right Givens scalars. FLA_Obj_create( dt_comp, min_m_n-1, n_GH, 0, 0, &G ); FLA_Obj_create( dt_comp, min_m_n-1, n_GH, 0, 0, &H ); // Create a real scaling factor. FLA_Obj_create( dt_real, 1, 1, 0, 0, &scale ); // Compute a scaling factor; If none is needed, sigma will be set to one. FLA_Svd_compute_scaling( A, scale ); // Scale the matrix if scale is non-unit. if ( !FLA_Obj_equals( scale, FLA_ONE ) ) FLA_Scal( scale, A ); if ( m_A < crossover_ratio * n_A ) { // Reduce the matrix to bidiagonal form. // Apply scalars to rotate elements on the superdiagonal to the real domain. // Extract the diagonal and superdiagonal from A. FLA_Bidiag_UT( A, T, S ); FLA_Bidiag_UT_realify( A, rL, rR ); FLA_Bidiag_UT_extract_real_diagonals( A, d, e ); // Form U and V. FLA_Bidiag_UT_form_U( A, T, U ); FLA_Bidiag_UT_form_V( A, S, V ); // Apply the realifying scalars in rL and rR to U and V, respectively. { FLA_Obj UL, UR; FLA_Obj VL, VR; FLA_Part_1x2( U, &UL, &UR, min_m_n, FLA_LEFT ); FLA_Part_1x2( V, &VL, &VR, min_m_n, FLA_LEFT ); FLA_Apply_diag_matrix( FLA_RIGHT, FLA_CONJUGATE, rL, UL ); FLA_Apply_diag_matrix( FLA_RIGHT, FLA_NO_CONJUGATE, rR, VL ); } // Perform a singular value decomposition on the bidiagonal matrix. r_val = FLA_Bsvd_v_opt_var1( n_iter_max, d, e, G, H, U, V, b_alg ); } else // if ( crossover_ratio * n_A <= m_A ) { FLA_Obj TQ, R; FLA_Obj AT, AB; FLA_Obj UL, UR; // Perform a QR factorization on A and form Q in U. FLA_QR_UT_create_T( A, &TQ ); FLA_QR_UT( A, TQ ); FLA_QR_UT_form_Q( A, TQ, U ); FLA_Obj_free( &TQ ); // Set the lower triangle of R to zero and then copy the upper // triangle of A to R. FLA_Part_2x1( A, &AT, &AB, n_A, FLA_TOP ); FLA_Obj_create( dt, n_A, n_A, 0, 0, &R ); FLA_Setr( FLA_LOWER_TRIANGULAR, FLA_ZERO, R ); FLA_Copyr( FLA_UPPER_TRIANGULAR, AT, R ); // Reduce the matrix to bidiagonal form. // Apply scalars to rotate elements on the superdiagonal to the real domain. // Extract the diagonal and superdiagonal from A. FLA_Bidiag_UT( R, T, S ); FLA_Bidiag_UT_realify( R, rL, rR ); FLA_Bidiag_UT_extract_real_diagonals( R, d, e ); // Form V from right Householder vectors in upper triangle of R. FLA_Bidiag_UT_form_V( R, S, V ); // Form U in R. FLA_Bidiag_UT_form_U( R, T, R ); // Apply the realifying scalars in rL and rR to U and V, respectively. FLA_Apply_diag_matrix( FLA_RIGHT, FLA_CONJUGATE, rL, R ); FLA_Apply_diag_matrix( FLA_RIGHT, FLA_NO_CONJUGATE, rR, V ); // Perform a singular value decomposition on the bidiagonal matrix. r_val = FLA_Bsvd_v_opt_var1( n_iter_max, d, e, G, H, R, V, b_alg ); // Multiply R into U, storing the result in A and then copying back // to U. FLA_Part_1x2( U, &UL, &UR, n_A, FLA_LEFT ); FLA_Gemm( FLA_NO_TRANSPOSE, FLA_NO_TRANSPOSE, FLA_ONE, UL, R, FLA_ZERO, A ); FLA_Copy( A, UL ); FLA_Obj_free( &R ); } // Copy the converged eigenvalues to the output vector. FLA_Copy( d, s ); // Sort the singular values and singular vectors in descending order. FLA_Sort_svd( FLA_BACKWARD, s, U, V ); // If the matrix was scaled, rescale the singular values. if ( !FLA_Obj_equals( scale, FLA_ONE ) ) FLA_Inv_scal( scale, s ); FLA_Obj_free( &scale ); FLA_Obj_free( &T ); FLA_Obj_free( &S ); FLA_Obj_free( &rL ); FLA_Obj_free( &rR ); FLA_Obj_free( &d ); FLA_Obj_free( &e ); FLA_Obj_free( &G ); FLA_Obj_free( &H ); return r_val; }
FLA_Error FLA_Hess_UT_step_unb_var2( FLA_Obj A, FLA_Obj T ) { FLA_Obj ATL, ATR, A00, a01, A02, ABL, ABR, a10t, alpha11, a12t, A20, a21, A22; FLA_Obj TTL, TTR, T00, t01, T02, TBL, TBR, t10t, tau11, t12t, T20, t21, T22; FLA_Obj yT, y0, yB, psi1, y2; FLA_Obj zT, z0, zB, zeta1, z2; FLA_Obj y, z; FLA_Obj inv_tau11; FLA_Obj minus_inv_tau11; FLA_Obj first_elem; FLA_Obj beta; FLA_Obj conj_beta; FLA_Obj dot_product; FLA_Obj a21_t, a21_b; FLA_Datatype datatype_A; dim_t m_A; dim_t b_alg; b_alg = FLA_Obj_length( T ); datatype_A = FLA_Obj_datatype( A ); m_A = FLA_Obj_length( A ); FLA_Obj_create( datatype_A, 1, 1, 0, 0, &inv_tau11 ); FLA_Obj_create( datatype_A, 1, 1, 0, 0, &minus_inv_tau11 ); FLA_Obj_create( datatype_A, 1, 1, 0, 0, &first_elem ); FLA_Obj_create( datatype_A, 1, 1, 0, 0, &beta ); FLA_Obj_create( datatype_A, 1, 1, 0, 0, &conj_beta ); FLA_Obj_create( datatype_A, 1, 1, 0, 0, &dot_product ); FLA_Obj_create( datatype_A, m_A, 1, 0, 0, &y ); FLA_Obj_create( datatype_A, m_A, 1, 0, 0, &z ); FLA_Part_2x2( A, &ATL, &ATR, &ABL, &ABR, 0, 0, FLA_TL ); FLA_Part_2x2( T, &TTL, &TTR, &TBL, &TBR, 0, 0, FLA_TL ); FLA_Part_2x1( y, &yT, &yB, 0, FLA_TOP ); FLA_Part_2x1( z, &zT, &zB, 0, FLA_TOP ); while ( FLA_Obj_length( ATL ) < b_alg ) { FLA_Repart_2x2_to_3x3( ATL, /**/ ATR, &A00, /**/ &a01, &A02, /* ************* */ /* ************************** */ &a10t, /**/ &alpha11, &a12t, ABL, /**/ ABR, &A20, /**/ &a21, &A22, 1, 1, FLA_BR ); FLA_Repart_2x2_to_3x3( TTL, /**/ TTR, &T00, /**/ &t01, &T02, /* ************* */ /* ************************** */ &t10t, /**/ &tau11, &t12t, TBL, /**/ TBR, &T20, /**/ &t21, &T22, 1, 1, FLA_BR ); FLA_Repart_2x1_to_3x1( yT, &y0, /* ** */ /* **** */ &psi1, yB, &y2, 1, FLA_BOTTOM ); FLA_Repart_2x1_to_3x1( zT, &z0, /* ** */ /* ***** */ &zeta1, zB, &z2, 1, FLA_BOTTOM ); /*------------------------------------------------------------*/ if ( FLA_Obj_length( A22 ) > 0 ) { FLA_Part_2x1( a21, &a21_t, &a21_b, 1, FLA_TOP ); // [ u21, tau11, a21 ] = House( a21 ); FLA_Househ2_UT( FLA_LEFT, a21_t, a21_b, tau11 ); // inv_tau11 = 1 / tau11; // minus_inv_tau11 = -1 / tau11; FLA_Set( FLA_ONE, inv_tau11 ); FLA_Inv_scalc( FLA_NO_CONJUGATE, tau11, inv_tau11 ); FLA_Copy( inv_tau11, minus_inv_tau11 ); FLA_Scal( FLA_MINUS_ONE, minus_inv_tau11 ); // Save first element of a21_t and set it to one so we can use a21 as // u21 in subsequent computations. We will restore a21_t later on. FLA_Copy( a21_t, first_elem ); FLA_Set( FLA_ONE, a21_t ); // y21 = A22' * u21; FLA_Gemv( FLA_CONJ_TRANSPOSE, FLA_ONE, A22, a21, FLA_ZERO, y2 ); // z21 = A22 * u21; FLA_Gemv( FLA_NO_TRANSPOSE, FLA_ONE, A22, a21, FLA_ZERO, z2 ); // beta = u21' * z21 / 2; // conj_beta = conj(beta); FLA_Dotc( FLA_CONJUGATE, a21, z2, beta ); FLA_Inv_scal( FLA_TWO, beta ); FLA_Copyt( FLA_CONJ_NO_TRANSPOSE, beta, conj_beta ); // y21' = ( y21' - beta / tau * u21' ) / tau; // y21 = ( y21 - conj(beta) / tau * u21 ) / tau; FLA_Scal( minus_inv_tau11, conj_beta ); FLA_Axpy( conj_beta, a21, y2 ); FLA_Scal( inv_tau11, y2 ); // z21 = ( z21 - beta / tau * u21 ) / tau; FLA_Scal( minus_inv_tau11, beta ); FLA_Axpy( beta, a21, z2 ); FLA_Scal( inv_tau11, z2 ); // a12t = a12t * ( I - u21 * u21' / tau ); // = a12t - ( a12t * u21 ) * u21' / tau; FLA_Dot( a12t, a21, dot_product ); FLA_Scal( minus_inv_tau11, dot_product ); FLA_Axpyt( FLA_CONJ_TRANSPOSE, dot_product, a21, a12t ); // A02 = A02 * ( I - u21 * u21' / tau ); // = A02 - ( A02 * u21 ) * u21' / tau; FLA_Gemv( FLA_NO_TRANSPOSE, FLA_ONE, A02, a21, FLA_ZERO, y0 ); FLA_Gerc( FLA_NO_CONJUGATE, FLA_CONJUGATE, minus_inv_tau11, y0, a21, A02 ); // A22 = A22 - u21 * y21' - z21 * u21'; FLA_Gerc( FLA_NO_CONJUGATE, FLA_CONJUGATE, FLA_MINUS_ONE, a21, y2, A22 ); FLA_Gerc( FLA_NO_CONJUGATE, FLA_CONJUGATE, FLA_MINUS_ONE, z2, a21, A22 ); // t01 = U20' * u21; FLA_Gemv( FLA_CONJ_TRANSPOSE, FLA_ONE, A20, a21, FLA_ZERO, t01 ); // Restore first element of a21. FLA_Copy( first_elem, a21_t ); } /*------------------------------------------------------------*/ FLA_Cont_with_3x3_to_2x2( &ATL, /**/ &ATR, A00, a01, /**/ A02, a10t, alpha11, /**/ a12t, /* ************** */ /* ************************ */ &ABL, /**/ &ABR, A20, a21, /**/ A22, FLA_TL ); FLA_Cont_with_3x3_to_2x2( &TTL, /**/ &TTR, T00, t01, /**/ T02, t10t, tau11, /**/ t12t, /* ************** */ /* ************************ */ &TBL, /**/ &TBR, T20, t21, /**/ T22, FLA_TL ); FLA_Cont_with_3x1_to_2x1( &yT, y0, psi1, /* ** */ /* **** */ &yB, y2, FLA_TOP ); FLA_Cont_with_3x1_to_2x1( &zT, z0, zeta1, /* ** */ /* ***** */ &zB, z2, FLA_TOP ); } FLA_Obj_free( &inv_tau11 ); FLA_Obj_free( &minus_inv_tau11 ); FLA_Obj_free( &first_elem ); FLA_Obj_free( &beta ); FLA_Obj_free( &conj_beta ); FLA_Obj_free( &dot_product ); FLA_Obj_free( &y ); FLA_Obj_free( &z ); return FLA_SUCCESS; }
FLA_Error FLA_Svd_compute_scaling( FLA_Obj A, FLA_Obj sigma ) { FLA_Datatype dt_real; FLA_Obj norm; FLA_Obj safmin; FLA_Obj prec; FLA_Obj rmin; FLA_Obj rmax; if ( FLA_Check_error_level() >= FLA_MIN_ERROR_CHECKING ) FLA_Svd_compute_scaling_check( A, sigma ); dt_real = FLA_Obj_datatype_proj_to_real( A ); FLA_Obj_create( dt_real, 1, 1, 0, 0, &norm ); FLA_Obj_create( dt_real, 1, 1, 0, 0, &prec ); FLA_Obj_create( dt_real, 1, 1, 0, 0, &safmin ); FLA_Obj_create( dt_real, 1, 1, 0, 0, &rmin ); FLA_Obj_create( dt_real, 1, 1, 0, 0, &rmax ); // Query safmin, precision. FLA_Mach_params( FLA_MACH_PREC, prec ); FLA_Mach_params( FLA_MACH_SFMIN, safmin ); //FLA_Obj_show( "safmin", safmin, "%20.12e", "" ); //FLA_Obj_show( "prec", prec, "%20.12e", "" ); // rmin = sqrt( safmin ) / prec; FLA_Copy( safmin, rmin ); FLA_Sqrt( rmin ); FLA_Inv_scal( prec, rmin ); // rmax = 1 / rmin; FLA_Copy( rmin, rmax ); FLA_Invert( FLA_NO_CONJUGATE, rmax ); //FLA_Obj_show( "rmin", rmin, "%20.12e", "" ); //FLA_Obj_show( "rmax", rmax, "%20.12e", "" ); // Find the maximum absolute value of A. FLA_Max_abs_value( A, norm ); if ( FLA_Obj_gt( norm, FLA_ZERO ) && FLA_Obj_lt( norm, rmin ) ) { // sigma = rmin / norm; FLA_Copy( rmin, sigma ); FLA_Inv_scal( norm, sigma ); } else if ( FLA_Obj_gt( norm, rmax ) ) { // sigma = rmax / norm; FLA_Copy( rmax, sigma ); FLA_Inv_scal( norm, sigma ); } else { // sigma = 1.0; FLA_Copy( FLA_ONE, sigma ); } FLA_Obj_free( &norm ); FLA_Obj_free( &prec ); FLA_Obj_free( &safmin ); FLA_Obj_free( &rmin ); FLA_Obj_free( &rmax ); return FLA_SUCCESS; }
FLA_Error FLA_Lyap_n_unb_var4( FLA_Obj isgn, FLA_Obj A, FLA_Obj C ) { FLA_Obj ATL, ATR, A00, a01, A02, ABL, ABR, a10t, alpha11, a12t, A20, a21, A22; FLA_Obj CTL, CTR, C00, c01, C02, CBL, CBR, c10t, gamma11, c12t, C20, c21, C22; FLA_Obj WTL, WTR, W00, w01, W02, WBL, WBR, w10t, omega11, w12t, W20, w21, W22; FLA_Obj W, omega; FLA_Scal( isgn, C ); FLA_Obj_create_conf_to( FLA_NO_TRANSPOSE, A, &W ); FLA_Obj_create( FLA_Obj_datatype( A ), 1, 1, 0, 0, &omega ); FLA_Part_2x2( A, &ATL, &ATR, &ABL, &ABR, 0, 0, FLA_BR ); FLA_Part_2x2( C, &CTL, &CTR, &CBL, &CBR, 0, 0, FLA_BR ); FLA_Part_2x2( W, &WTL, &WTR, &WBL, &WBR, 0, 0, FLA_BR ); while ( FLA_Obj_length( CTL ) > 0 ){ FLA_Repart_2x2_to_3x3( ATL, /**/ ATR, &A00, &a01, /**/ &A02, &a10t, &alpha11, /**/ &a12t, /* ************* */ /* ************************** */ ABL, /**/ ABR, &A20, &a21, /**/ &A22, 1, 1, FLA_TL ); FLA_Repart_2x2_to_3x3( CTL, /**/ CTR, &C00, &c01, /**/ &C02, &c10t, &gamma11, /**/ &c12t, /* ************* */ /* ************************** */ CBL, /**/ CBR, &C20, &c21, /**/ &C22, 1, 1, FLA_TL ); FLA_Repart_2x2_to_3x3( WTL, /**/ WTR, &W00, &w01, /**/ &W02, &w10t, &omega11, /**/ &w12t, /* ************* */ /* ************************** */ WBL, /**/ WBR, &W20, &w21, /**/ &W22, 1, 1, FLA_TL ); /*------------------------------------------------------------*/ // gamma11 = gamma11 / ( alpha11 + alpha11' ); FLA_Copyt( FLA_CONJ_NO_TRANSPOSE, alpha11, omega ); FLA_Mult_add( FLA_ONE, alpha11, omega ); FLA_Inv_scal( omega, gamma11 ); // c01 = c01 - a01 * gamma11; FLA_Axpys( FLA_MINUS_ONE, gamma11, a01, FLA_ONE, c01 ); // c01 = inv( triu(A00) + conj(alpha) * I ) * c01; FLA_Copyrt( FLA_UPPER_TRIANGULAR, FLA_NO_TRANSPOSE, A00, W00 ); FLA_Shift_diag( FLA_CONJUGATE, alpha11, W00 ); FLA_Trsv( FLA_UPPER_TRIANGULAR, FLA_NO_TRANSPOSE, FLA_NONUNIT_DIAG, W00, c01 ); // C00 = C00 - a01 * c01' - c01 * a01'; FLA_Her2( FLA_UPPER_TRIANGULAR, FLA_MINUS_ONE, a01, c01, C00 ); /*------------------------------------------------------------*/ FLA_Cont_with_3x3_to_2x2( &ATL, /**/ &ATR, A00, /**/ a01, A02, /* ************** */ /* ************************ */ a10t, /**/ alpha11, a12t, &ABL, /**/ &ABR, A20, /**/ a21, A22, FLA_BR ); FLA_Cont_with_3x3_to_2x2( &CTL, /**/ &CTR, C00, /**/ c01, C02, /* ************** */ /* ************************ */ c10t, /**/ gamma11, c12t, &CBL, /**/ &CBR, C20, /**/ c21, C22, FLA_BR ); FLA_Cont_with_3x3_to_2x2( &WTL, /**/ &WTR, W00, /**/ w01, W02, /* ************** */ /* ************************ */ w10t, /**/ omega11, w12t, &WBL, /**/ &WBR, W20, /**/ w21, W22, FLA_BR ); } FLA_Obj_free( &W ); FLA_Obj_free( &omega ); return FLA_SUCCESS; }
int main( int argc, char** argv ) { FLA_Datatype comptype = COMPTYPE; FLA_Datatype realtype = REALTYPE; dim_t m; FLA_Obj a, aT, aB, a0, a1, a2; FLA_Obj v, vT, vB, v0, v1, v2; FLA_Error init_result; int use_abs = 1; if ( argc == 3 ) { m = atoi(argv[1]); use_abs = atoi(argv[2]); } else { fprintf(stderr, " \n"); fprintf(stderr, "Usage: %s m use_abs\n", argv[0]); fprintf(stderr, " m : test vector length\n"); fprintf(stderr, " use_abs : 0 - norm (realtype), 1 - abs (complex type)\n"); fprintf(stderr, " \n"); return -1; } if ( m == 0 ) return 0; FLA_Init_safe( &init_result ); FLA_Obj_create( comptype, m, 1, 0, 0, &a ); FLA_Obj_create( use_abs ? comptype : realtype, m, 1, 0, 0, &v ); FLA_Random_matrix( a ); FLA_Set( FLA_ZERO, v ); FLA_Obj_fshow( stdout, "- a -", a, "% 6.4e", "--" ); // Normalize a vector FLA_Part_2x1( a, &aT, &aB, 0, FLA_TOP ); FLA_Part_2x1( v, &vT, &vB, 0, FLA_TOP ); while ( FLA_Obj_length( aB ) > 0 ) { FLA_Repart_2x1_to_3x1( aT, &a0, &a1, aB, &a2, 1, FLA_BOTTOM ); FLA_Repart_2x1_to_3x1( vT, &v0, &v1, vB, &v2, 1, FLA_BOTTOM ); // -------------------------------------------- if ( use_abs ) { // a and v are complex datatype FLA_Copy( a1, v1 ); FLA_Absolute_value( v1 ); } else { // v is real datatype FLA_Nrm2( a1, v1 ); } if ( FLA_Obj_equals( v1, FLA_ZERO ) ) printf( " ZERO DETECTED\n" ); else FLA_Inv_scal( v1, a1 ); // Normalize the scalar // -------------------------------------------- FLA_Cont_with_3x1_to_2x1( &aT, a0, a1, &aB, a2, FLA_TOP ); FLA_Cont_with_3x1_to_2x1( &vT, v0, v1, &vB, v2, FLA_TOP ); } FLA_Obj_fshow( stdout, "- a -", a, "% 6.4e", "--" ); FLA_Obj_fshow( stdout, "- v -", v, "% 6.4e", "--" ); // Check whether it is normalized FLA_Part_2x1( a, &aT, &aB, 0, FLA_TOP ); FLA_Part_2x1( v, &vT, &vB, 0, FLA_TOP ); while ( FLA_Obj_length( aB ) > 0 ) { FLA_Repart_2x1_to_3x1( aT, &a0, &a1, aB, &a2, 1, FLA_BOTTOM ); FLA_Repart_2x1_to_3x1( vT, &v0, &v1, vB, &v2, 1, FLA_BOTTOM ); // -------------------------------------------- if ( use_abs ) { // a and v are same datatype FLA_Copy( a1, v1 ); FLA_Absolute_value( v1 ); } else { // v is realdatatype FLA_Nrm2( a1, v1 ); } // -------------------------------------------- FLA_Cont_with_3x1_to_2x1( &aT, a0, a1, &aB, a2, FLA_TOP ); FLA_Cont_with_3x1_to_2x1( &vT, v0, v1, &vB, v2, FLA_TOP ); } FLA_Obj_fshow( stdout, " - all should be one - ", v, "% 6.4e", "--"); FLA_Obj_free( &a ); FLA_Obj_free( &v ); FLA_Finalize_safe( init_result ); }
FLA_Error FLA_Fill_with_logarithmic_dist( FLA_Obj alpha, FLA_Obj x ) { FLA_Obj lT, l0, lB, lambda1, l2; FLA_Obj l, k, alpha2; FLA_Datatype dt_real; dim_t n_x; if ( FLA_Check_error_level() >= FLA_MIN_ERROR_CHECKING ) FLA_Fill_with_logarithmic_dist_check( alpha, x ); dt_real = FLA_Obj_datatype_proj_to_real( x ); n_x = FLA_Obj_vector_dim( x ); // Create a local counter to increment as we create the distribution. FLA_Obj_create( dt_real, 1, 1, 0, 0, &k ); // Create a local vector l. We will work with this vector, which is // the same length as x, so that we can use vertical partitioning. FLA_Obj_create( dt_real, n_x, 1, 0, 0, &l ); // Create a local real scalar alpha2 of the same precision as // alpha. Then copy alpha to alpha2, which will convert the // complex value to real, if necessary (ie: if alpha is complex). FLA_Obj_create( dt_real, 1, 1, 0, 0, &alpha2 ); FLA_Copy( alpha, alpha2 ); // Initialize k to 0. FLA_Set( FLA_ZERO, k ); FLA_Part_2x1( l, &lT, &lB, 0, FLA_TOP ); while ( FLA_Obj_length( lB ) > 0 ) { FLA_Repart_2x1_to_3x1( lT, &l0, /* ** */ /* ******* */ &lambda1, lB, &l2, 1, FLA_BOTTOM ); /*------------------------------------------------------------*/ // lambda1 = alpha^k; FLA_Pow( alpha2, k, lambda1 ); // k = k + 1; FLA_Mult_add( FLA_ONE, FLA_ONE, k ); /*------------------------------------------------------------*/ FLA_Cont_with_3x1_to_2x1( &lT, l0, lambda1, /* ** */ /* ******* */ &lB, l2, FLA_TOP ); } // Normalize by last element. FLA_Part_2x1( l, &lT, &lB, 1, FLA_BOTTOM ); FLA_Inv_scal( lB, l ); // Overwrite x with the distribution we created in l. FLA_Copy( l, x ); FLA_Obj_free( &l ); FLA_Obj_free( &k ); FLA_Obj_free( &alpha2 ); return FLA_SUCCESS; }