@@ -576,6 +576,299 @@ def expm1(x, prec)
576576 exp_prec > 0 ? exp ( x , exp_prec ) . sub ( 1 , prec ) : BigDecimal ( -1 )
577577 end
578578
579+ # call-seq:
580+ # erf(decimal, numeric) -> BigDecimal
581+ #
582+ # Computes the error function of +decimal+ to the specified number of digits of
583+ # precision, +numeric+.
584+ #
585+ # If +decimal+ is NaN, returns NaN.
586+ #
587+ # BigMath.erf(BigDecimal('1'), 32).to_s
588+ # #=> "0.84270079294971486934122063508261e0"
589+ #
590+ def erf ( x , prec )
591+ prec = BigDecimal ::Internal . coerce_validate_prec ( prec , :erf )
592+ x = BigDecimal ::Internal . coerce_to_bigdecimal ( x , prec , :erf )
593+ return BigDecimal ::Internal . nan_computation_result if x . nan?
594+ return BigDecimal ( x . infinite? ) if x . infinite?
595+ return BigDecimal ( 0 ) if x == 0
596+ return -erf ( -x , prec ) if x < 0
597+ return BigDecimal ( 1 ) if x > 5000000000 # erf(5000000000) > 1 - 1e-10000000000000000000
598+
599+ if x > 8
600+ xf = x . to_f
601+ log10_erfc = -xf ** 2 / Math . log ( 10 ) - Math . log10 ( xf * Math ::PI ** 0.5 )
602+ erfc_prec = [ prec + log10_erfc . ceil , 1 ] . max
603+ erfc = _erfc_bit_burst ( x , erfc_prec + BigDecimal . double_fig )
604+ return BigDecimal ( 1 ) . sub ( erfc , prec ) if erfc
605+ end
606+
607+ _erf_bit_burst ( x , prec + BigDecimal . double_fig ) . mult ( 1 , prec )
608+ end
609+
610+ def erfc ( x , prec )
611+ prec = BigDecimal ::Internal . coerce_validate_prec ( prec , :erfc )
612+ x = BigDecimal ::Internal . coerce_to_bigdecimal ( x , prec , :erfc )
613+ return BigDecimal ::Internal . nan_computation_result if x . nan?
614+ return BigDecimal ( 1 - x . infinite? ) if x . infinite?
615+ return BigDecimal ( 1 ) . sub ( erf ( x , prec ) , prec ) if x < 0
616+ return BigDecimal ( 0 ) if x > 5000000000 # erfc(5000000000) < 1e-10000000000000000000 (underflow)
617+
618+ if x >= 8
619+ y = _erfc_bit_burst ( x , prec + BigDecimal . double_fig )
620+ return y . mult ( 1 , prec ) if y
621+ end
622+
623+ # erfc(x) = 1 - erf(x) < exp(-x**2)/x/sqrt(pi)
624+ # Precision of erf(x) needs about log10(exp(-x**2)) extra digits
625+ log10 = 2.302585092994046
626+ high_prec = prec + BigDecimal . double_fig + ( x . to_f **2 / log10 ) . ceil
627+ BigDecimal ( 1 ) . sub ( _erf_bit_burst ( x , high_prec ) , prec )
628+ end
629+
630+ # Calculates erf(x) using bit-burst algorithm.
631+ private_class_method def _erf_bit_burst ( x , prec )
632+ x = BigDecimal ::Internal . coerce_to_bigdecimal ( x , prec , :erf )
633+ prec = BigDecimal ::Internal . coerce_validate_prec ( prec , :erf )
634+
635+ return BigDecimal ( 0 ) if x > 5000000000 # erfc underflows
636+ x = x . mult ( 1 , [ prec - ( x . to_f **2 /Math . log ( 10 ) ) . floor , 1 ] . max )
637+
638+ calculated_x = BigDecimal ( 0 )
639+ erf_exp2 = BigDecimal ( 0 )
640+ digits = 8
641+ xf = x . to_f
642+ scale = 2 * exp ( -x . mult ( x , prec ) , prec ) . div ( PI ( prec ) . sqrt ( prec ) , prec )
643+
644+ until x . zero?
645+ partial = x . truncate ( digits )
646+ digits *= 2
647+ next if partial . zero?
648+
649+ erf_exp2 = _erf_exp2_binary_splitting ( partial , calculated_x , erf_exp2 , prec )
650+ calculated_x += partial
651+ x -= partial
652+ end
653+ erf_exp2 . mult ( scale , prec )
654+ end
655+
656+ # Calculates erfc(x) using bit-burst algorithm.
657+ private_class_method def _erfc_bit_burst ( x , prec ) # :nodoc:
658+ digits = ( x . exponent + 1 ) * 40
659+
660+ calculated_x = x . truncate ( digits )
661+ f = _erfc_exp2_asymptotic_binary_splitting ( calculated_x , prec )
662+ return unless f
663+
664+ scale = 2 * exp ( -x . mult ( x , prec ) , prec ) . div ( PI ( prec ) . sqrt ( prec ) , prec )
665+ x -= calculated_x
666+
667+ until x . zero?
668+ digits *= 2
669+ partial = x . truncate ( digits )
670+ next if partial . zero?
671+
672+ f = _erfc_exp2_inv_inv_binary_splitting ( partial , calculated_x , f , prec )
673+ calculated_x += partial
674+ x -= partial
675+ end
676+ f . mult ( scale , prec )
677+ end
678+
679+ # Matrix multiplication for binary splitting method in erf/erfc calculation
680+ private_class_method def _bs_matrix_mult ( m1 , m2 , size , prec ) # :nodoc:
681+ ( size * size ) . times . map do |i |
682+ size . times . map do |k |
683+ m1 [ i / size * size + k ] . mult ( m2 [ size * k + i % size ] , prec )
684+ end . reduce { |a , b | a . add ( b , prec ) }
685+ end
686+ end
687+
688+ # Matrix/Vector weighted sum for binary splitting method in erf/erfc calculation
689+ private_class_method def _bs_weighted_sum ( m1 , w1 , m2 , w2 , prec ) # :nodoc:
690+ m1 . zip ( m2 ) . map { |v1 , v2 | ( v1 * w1 ) . add ( v2 * w2 , prec ) }
691+ end
692+
693+ # Calculates Taylor expansion of erf(x+a)*exp((x+a)**2)*sqrt(pi)/2 with binary splitting method.
694+ private_class_method def _erf_exp2_binary_splitting ( x , a , f_a , prec ) # :nodoc:
695+ cexponent = Math . log10 ( [ 2 * a , Math . sqrt ( 2 ) ] . max . to_f ) + log10 ( x . abs , 10 )
696+ log10f = Math . log ( 10 )
697+
698+ steps = BigDecimal . save_exception_mode do
699+ BigDecimal . mode ( BigDecimal ::EXCEPTION_UNDERFLOW , false )
700+ ( 2 ..) . bsearch do |n |
701+ x . to_f ** 2 < n && n * cexponent + Math . lgamma ( n / 2 ) [ 0 ] / log10f + n * Math . log10 ( 2 ) - Math . lgamma ( n - 1 ) [ 0 ] / log10f < -prec + x . to_f **2 / log10f
702+ end
703+ end
704+
705+ if a == 0
706+ return x . mult ( 1 + BigDecimal ::Internal . taylor_sum_binary_splitting ( 2 * x * x , ( steps / 2 ) . times . map { 2 * _1 + 3 } , prec ) , prec )
707+ end
708+
709+ # First, calculate a matrix that represents the sum of the Taylor series:
710+ # SumMatrix = (((((...+I)x*M4+I)*x*M3+I)*M2*x+I)*M1*x+I)
711+ # Where Mi is a 2x2 matrix that generates the next coefficients of Taylor series:
712+ # Vector(c4, c5) = M4*M3*M2*M1*Vector(c0, c1)
713+ # And then calculates:
714+ # SumMatrix * Vector(c0, c1) = Vector(c0+c1*x+c2*x**2+..., _)
715+ # In this binary splitting method, adjacent two operations are combined into one repeatedly.
716+ # ((...) * x * A + B) / C is the form of each operation. A and B are 2x2 matrices, C is a scalar.
717+ zero = BigDecimal ( 0 )
718+ two = BigDecimal ( 2 )
719+ two_a = two * a
720+ operations = steps . times . map do |i |
721+ n = BigDecimal ( 2 + i )
722+ [ [ zero , n , two , two_a ] , [ n , zero , zero , n ] , n ]
723+ end
724+
725+ while operations . size > 1
726+ xpow = xpow ? xpow . mult ( xpow , prec ) : x . mult ( 1 , prec )
727+ operations = operations . each_slice ( 2 ) . map do |op1 , op2 |
728+ # Combine two operations into one:
729+ # (((Remaining * x * A2 + B2) / C2) * x * A1 + B1) / C1
730+ # ((Remaining * (x*x) * (A2*A1) + (x*B2*A1+B1*C2)) / (C1*C2)
731+ # Therefore, combined operation can be represented as:
732+ # Anext = A2 * A1
733+ # Bnext = x * B2 * A1 + B1 * C2
734+ # Cnext = C1 * C2
735+ # xnext = x * x
736+ a1 , b1 , c1 = op1
737+ a2 , b2 , c2 = op2 || [ [ zero ] * 4 , [ zero ] * 4 , BigDecimal ( 1 ) ]
738+ [
739+ _bs_matrix_mult ( a2 , a1 , 2 , prec ) ,
740+ _bs_weighted_sum ( _bs_matrix_mult ( b2 , a1 , 2 , prec ) , xpow , b1 , c2 , prec ) ,
741+ c1 . mult ( c2 , prec ) ,
742+ ]
743+ end
744+ end
745+ _ , sum_matrix , denominator = operations . first
746+ ( sum_matrix [ 1 ] + f_a * ( 2 * a * sum_matrix [ 1 ] + sum_matrix [ 0 ] ) ) . div ( denominator , prec )
747+ end
748+
749+ # Calculates asymptotic expansion of erfc(x)*exp(x**2)*sqrt(pi)/2 with binary splitting method
750+ private_class_method def _erfc_exp2_asymptotic_binary_splitting ( x , prec ) # :nodoc:
751+ # Let f(x) = erfc(x)*sqrt(pi)*exp(x**2)/2
752+ # f(x) satisfies the following differential equation:
753+ # 2*x*f(x) = f'(x) + 1
754+ # From the above equation, we can derive the following asymptotic expansion:
755+ # f(x) = (0..kmax).sum { (-1)**k * (2*k)! / 4**k / k! / x**(2*k)) } / x
756+
757+ # This asymptotic expansion does not converge.
758+ # But if there is a k that satisfies (2*k)! / 4**k / k! / x**(2*k) < 10**(-prec),
759+ # It is enough to calculate erfc within the given precision.
760+ # Using Stirling's approximation, we can simplify this condition to:
761+ # sqrt(2)/2 + k*log(k) - k - 2*k*log(x) < -prec*log(10)
762+ # and the left side is minimized when k = x**2.
763+ xf = x . to_f
764+ kmax = ( 1 ..( xf ** 2 ) . floor ) . bsearch do |k |
765+ Math . log ( 2 ) / 2 + k * Math . log ( k ) - k - 2 * k * Math . log ( xf ) < -prec * Math . log ( 10 )
766+ end
767+ return unless kmax
768+
769+ # Convert asymptotic expansion to nested form:
770+ # 1 + a/x + a*b/x/x + a*b*c/x/x/x + a*b*c/x/x/x*rest
771+ # = 1 + (a/x) * (1 + (b/x) * (1 + (c/x) * (1 + rest)))
772+ #
773+ # And calculate it with binary splitting:
774+ # (a1/d + b1/d * (a2/d + b2/d * (rest)))
775+ # = ((a1*d+b1*a2)/(d*d) + b1*b2/(d*denominator) * (rest)))
776+ denominator = x . mult ( x , prec ) . mult ( 2 , prec )
777+ fractions = ( 1 ..kmax ) . map do |k |
778+ [ denominator , BigDecimal ( 1 - 2 * k ) ]
779+ end
780+ while fractions . size > 1
781+ fractions = fractions . each_slice ( 2 ) . map do |fraction1 , fraction2 |
782+ a1 , b1 = fraction1
783+ a2 , b2 = fraction2 || [ BigDecimal ( 0 ) , denominator ]
784+ [
785+ a1 . mult ( denominator , prec ) . add ( b1 . mult ( a2 , prec ) , prec ) ,
786+ b1 . mult ( b2 , prec ) ,
787+ ]
788+ end
789+ denominator = denominator . mult ( denominator , prec )
790+ end
791+ sum = fractions [ 0 ] [ 0 ] . add ( fractions [ 0 ] [ 1 ] , prec ) . div ( denominator , prec )
792+ sum . div ( x , prec ) / 2
793+ end
794+
795+ # Calculates f(1/(a + x)) where f(x) = (sqrt(pi)/2) * exp(1/x**2) * erfc(1/x)
796+ # f(1/(a+x)) = f(1/a - x/(a*(a+x)))
797+ private_class_method def _erfc_exp2_inv_inv_binary_splitting ( x , a , f_inva , prec )
798+ return f_inva if x . zero?
799+ # f(x) satisfies the following differential equation:
800+ # (1/a+w)**3*f'(1/a+w) + 2*f(1/a+w) = 1/a + w
801+ # From the above equation, we can derive the following Taylor expansion around x=a:
802+ # Coefficients: f(1/a + w) = c0 + c1*w + c2*w**2 + c3*w**3 + ...
803+ # Constraints:
804+ # (w**3 + 3*w**2/a + 3*w/a**2 + 1/a**3) * (c1 + 2*c2*w + 3*c3*w**2 + 4*c4*w**3 + ...)
805+ # + 2 * (c0 + c1*w + c2*w**2 + c3*w**3 + ...) = 1/a + w
806+ # Recurrence relations:
807+ # c0 = f(1/a)
808+ # c1 = a**2 - 2*c0*a**3
809+ # c2 = (a**3 - 3*c1*a - 2*c1*a**3) / 2
810+ # c3 = -(3*c1*a**2 + 6*c2*a + 2*c2*a**3) / 3
811+ # c(n) = -((n-3)*c(n-3)*a**3 + 3*(n-2)*c(n-2)*a**2 + 3*(n-1)*c(n-1)*a + 2*c(n-1)*a**3) / n
812+
813+ aa = a . mult ( a , prec )
814+ aaa = aa . mult ( a , prec )
815+ c0 = f_inva
816+ c1 = ( aa - 2 * c0 * aaa ) . mult ( 1 , prec )
817+ c2 = ( aaa - 3 * c1 * a - 2 * c1 * aaa ) . div ( 2 , prec )
818+
819+ # Estimate the number of steps needed to achieve the required precision
820+ low_prec = 16
821+ w = x . div ( a . mult ( a + x , low_prec ) , low_prec )
822+ wpow = w . mult ( w , low_prec )
823+ cm3 , cm2 , cm1 = [ c0 , c1 , c2 ] . map { _1 . mult ( 1 , low_prec ) }
824+ a_low , aa_low , aaa_low = [ a , aa , aaa ] . map { _1 . mult ( 1 , low_prec ) }
825+ step = ( 3 ..) . find do |n |
826+ wpow = wpow . mult ( w , low_prec )
827+ cn = -( ( n - 3 ) * cm3 * aaa_low + 3 * aa_low * ( n - 2 ) * cm2 + 3 * a_low * ( n - 1 ) * cm1 + 2 * cm1 * aaa_low ) . div ( n , low_prec )
828+ cm3 , cm2 , cm1 = cm2 , cm1 , cn
829+ cn . mult ( wpow , low_prec ) . exponent < -prec
830+ end
831+
832+ # Let M(n) be a 3x3 matrix that transforms (c(n),c(n+1),c(n+2)) to (c(n-1),c(n),c(n+1))
833+ # Mn = | 0 1 0 |
834+ # | 0 0 1 |
835+ # | -(n-3)*aaa/n -3*(n-2)*aa/n -2*aaa-3*(n-1)*a/n |
836+ # Vector(c6,c7,c8) = M6*M5*M4*M3*M2*M1 * Vector(c0,c1,c2)
837+ # Vector(c0+c1*y/z+c2*(y/z)**2+..., _, _) = (((... + I)*M3*y/z + I)*M2*y/z + I)*M1*y/z + I) * Vector(c2, c1, c0)
838+ # Perform binary splitting on this nested parenthesized calculation by using the following formula:
839+ # (((...)*A2*y/z + B2)/D2 * A1*y/z + B1)/D1 = (((...)*(A2*A1)*(y*y)/z + (B2*A1*y+z*D2*B1)) / (D1*D2*z)
840+ # where A_n, Bn are matrices and Dn are scalars
841+
842+ zero = BigDecimal ( 0 )
843+ operations = ( 3 ..step + 2 ) . map do |n |
844+ bign = BigDecimal ( n )
845+ [
846+ [
847+ zero , bign , zero ,
848+ zero , zero , bign ,
849+ BigDecimal ( -( n - 3 ) * aaa ) , -3 * ( n - 2 ) * aa , -2 * aaa - 3 * ( n - 1 ) * a
850+ ] ,
851+ [ bign , zero , zero , zero , bign , zero , zero , zero , bign ] ,
852+ bign
853+ ]
854+ end
855+
856+ z = a . mult ( a + x , prec )
857+ while operations . size > 1
858+ y = y ? y . mult ( y , prec ) : -x . mult ( 1 , prec )
859+ operations = operations . each_slice ( 2 ) . map do |op1 , op2 |
860+ a1 , b1 , d1 = op1
861+ a2 , b2 , d2 = op2 || [ [ zero ] * 9 , [ zero ] * 9 , BigDecimal ( 1 ) ]
862+ [
863+ _bs_matrix_mult ( a2 , a1 , 3 , prec ) ,
864+ _bs_weighted_sum ( _bs_matrix_mult ( b2 , a1 , 3 , prec ) , y , b1 , d2 . mult ( z , prec ) , prec ) ,
865+ d1 . mult ( d2 , prec ) . mult ( z , prec ) ,
866+ ]
867+ end
868+ end
869+ _ , sum_matrix , denominator = operations [ 0 ]
870+ ( sum_matrix [ 0 ] * c0 + sum_matrix [ 1 ] * c1 + sum_matrix [ 2 ] * c2 ) . div ( denominator , prec )
871+ end
579872
580873 # call-seq:
581874 # PI(numeric) -> BigDecimal
@@ -588,38 +881,18 @@ def expm1(x, prec)
588881 #
589882 def PI ( prec )
590883 prec = BigDecimal ::Internal . coerce_validate_prec ( prec , :PI )
591- n = prec + BigDecimal . double_fig
592- zero = BigDecimal ( "0" )
593- one = BigDecimal ( "1" )
594- two = BigDecimal ( "2" )
595-
596- m25 = BigDecimal ( "-0.04" )
597- m57121 = BigDecimal ( "-57121" )
598-
599- pi = zero
600-
601- d = one
602- k = one
603- t = BigDecimal ( "-80" )
604- while d . nonzero? && ( ( m = n - ( pi . exponent - d . exponent ) . abs ) > 0 )
605- m = BigDecimal . double_fig if m < BigDecimal . double_fig
606- t = t *m25
607- d = t . div ( k , m )
608- k = k +two
609- pi = pi + d
610- end
611-
612- d = one
613- k = one
614- t = BigDecimal ( "956" )
615- while d . nonzero? && ( ( m = n - ( pi . exponent - d . exponent ) . abs ) > 0 )
616- m = BigDecimal . double_fig if m < BigDecimal . double_fig
617- t = t . div ( m57121 , n )
618- d = t . div ( k , m )
619- pi = pi + d
620- k = k +two
884+ n = prec + BigDecimal . double_fig
885+ a = BigDecimal ( 1 )
886+ b = BigDecimal ( 0.5 , 0 ) . sqrt ( n )
887+ s = BigDecimal ( 0.25 , 0 )
888+ t = 1
889+ while a != b && ( a - b ) . exponent > 1 - n
890+ c = ( a - b ) . div ( 2 , n )
891+ a , b = ( a + b ) . div ( 2 , n ) , ( a * b ) . sqrt ( n )
892+ s = s . sub ( c * c * t , n )
893+ t *= 2
621894 end
622- pi . mult ( 1 , prec )
895+ ( a * b ) . div ( s , prec )
623896 end
624897
625898 # call-seq:
0 commit comments