def _evalf_(self, x, parent=None, algorithm=None): """ EXAMPLES:: sage: arcsec(2).n(100) 1.0471975511965977461542144611 sage: arcsec(1/2).n(100) NaN TESTS: Test complex input:: sage: arcsec(complex(1,1)) # rel tol 1e-15 (1.118517879643706+0.5306375309525178j) """ if parent is float: return math.acos(1 / x) try: return (1 / x).arccos() except AttributeError: # Usually this means that x is of type 'complex' from sage.rings.complex_double import CDF return complex(CDF(1 / x).arccos())
def _evalf_(self, x, parent=None, algorithm=None): """ EXAMPLES:: sage: arccsc(2).n(100) 0.52359877559829887307710723055 sage: float(arccsc(2)) 0.52359877559829... TESTS: Test complex input:: sage: arccsc(complex(1,1)) # rel tol 1e-15 (0.45227844715119064-0.5306375309525178j) """ if parent is float: return math.asin(1 / x) try: return (1 / x).arcsin() except AttributeError: # Usually this means that x is of type 'complex' from sage.rings.complex_double import CDF return complex(CDF(1 / x).arcsin())
def real(x): """ Return the real part of x. EXAMPLES:: sage: z = 1+2*I sage: real(z) 1 sage: real(5/3) 5/3 sage: a = 2.5 sage: real(a) 2.50000000000000 sage: type(real(a)) <type 'sage.rings.real_mpfr.RealLiteral'> """ #Try to all the .real() method try: return x.real() except AttributeError: pass #Try to coerce x into RDF. If that #succeeds, then we can just return x try: rdf_x = RDF(x) return x except TypeError: pass #Finall try to coerce x into CDF and call #the .real() method. return CDF(x).real()
def _evalf_(self, x, parent=None, algorithm=None): """ EXAMPLES:: sage: arccot(1/2).n(100) 1.1071487177940905030170654602 sage: float(arccot(1/2)) 1.1071487177940904 TESTS: Test complex input:: sage: arccot(complex(1,1)) # rel tol 1e-15 (0.5535743588970452-0.4023594781085251j) """ if parent is float: return math.pi / 2 - math.atan(x) from sage.symbolic.constants import pi try: return parent(pi / 2 - x.arctan()) except AttributeError: # Usually this means that x is of type 'complex' from sage.rings.complex_double import CDF return complex(pi / 2 - CDF(x).arctan())
def sqrt(x): """ Returns a square root of x. This function (``numerical_sqrt``) is deprecated. Use ``sqrt(x, prec=n)`` instead. EXAMPLES:: sage: numerical_sqrt(10.1) doctest:1: DeprecationWarning: numerical_sqrt is deprecated, use sqrt(x, prec=n) instead See http://trac.sagemath.org/5404 for details. 3.17804971641414 sage: numerical_sqrt(9) 3 """ from sage.misc.superseded import deprecation deprecation(5404, "numerical_sqrt is deprecated, use sqrt(x, prec=n) instead") try: return x.sqrt() except (AttributeError, ValueError): try: return RDF(x).sqrt() except TypeError: return CDF(x).sqrt()
def imag(x): """ Return the imaginary part of x. """ try: return x.imag() except AttributeError: return CDF(x).imag()
def arg(x): """ Return the argument of a complex number `x`. EXAMPLES:: sage: z = CC(1,2) sage: theta = arg(z) sage: cos(theta)*abs(z) 1.00000000000000 sage: sin(theta)*abs(z) 2.00000000000000 """ try: return x.arg() except AttributeError: return CDF(x).arg()
def sqrt(x): """ Return a square root of x. EXAMPLES:: sage: sqrt(10.1) 3.17804971641414 sage: sqrt(9) 3 """ try: return x.sqrt() except (AttributeError, ValueError): try: return RDF(x).sqrt() except TypeError: return CDF(x).sqrt()
def eta(x): r""" Returns the value of the eta function at `x`, which must be in the upper half plane. The `\eta` function is .. MATH:: \eta(z) = e^{\pi i z / 12} \prod_{n=1}^{\infty}(1-e^{2\pi inz}) EXAMPLES:: sage: eta(1+I) 0.7420487758365647 + 0.1988313702299107*I """ try: return x.eta() except AttributeError: return CDF(x).eta()
def _evalf_(self, n, z, parent=None, algorithm=None): """ EXAMPLES:: sage: N(lambert_w(1)) 0.567143290409784 sage: lambert_w(RealField(100)(1)) 0.56714329040978387299996866221 SciPy is used to evaluate for float, RDF, and CDF inputs:: sage: lambert_w(RDF(1)) 0.5671432904097838 sage: lambert_w(float(1)) 0.5671432904097838 sage: lambert_w(CDF(1)) 0.5671432904097838 sage: lambert_w(complex(1)) (0.5671432904097838+0j) sage: lambert_w(RDF(-1)) # abs tol 2e-16 -0.31813150520476413 + 1.3372357014306895*I sage: lambert_w(float(-1)) # abs tol 2e-16 (-0.31813150520476413+1.3372357014306895j) """ R = parent or s_parent(z) if R is float or R is RDF: from scipy.special import lambertw res = lambertw(z, n) # SciPy always returns a complex value, make it real if possible if not res.imag: return R(res.real) elif R is float: return complex(res) else: return CDF(res) elif R is complex or R is CDF: from scipy.special import lambertw return R(lambertw(z, n)) else: import mpmath return mpmath_utils.call(mpmath.lambertw, z, n, parent=R)
def sqrt(x): """ Returns a square root of x. This function (``numerical_sqrt``) is deprecated. Use ``sqrt(x, prec=n)`` instead. EXAMPLES:: sage: numerical_sqrt(10.1) doctest:1: DeprecationWarning: numerical_sqrt is deprecated, use sqrt(x, prec=n) instead 3.17804971641414 sage: numerical_sqrt(9) 3 """ from sage.misc.misc import deprecation deprecation("numerical_sqrt is deprecated, use sqrt(x, prec=n) instead") try: return x.sqrt() except (AttributeError, ValueError): try: return RDF(x).sqrt() except TypeError: return CDF(x).sqrt()
def __init__(self, x, universe=None, check=True, immutable=False, cr=False, cr_str=None, use_sage_types=False): """ Create a sequence. EXAMPLES:: sage: Sequence([1..5]) [1, 2, 3, 4, 5] sage: a = Sequence([1..3], universe=QQ, check=False, immutable=True, cr=True, cr_str=False, use_sage_types=True) sage: a [ 1, 2, 3 ] sage: a = Sequence([1..5], universe=QQ, check=False, immutable=True, cr_str=True, use_sage_types=True) sage: a [1, 2, 3, 4, 5] sage: a._Sequence__cr_str True sage: a.__str__() '[\n1,\n2,\n3,\n4,\n5\n]' """ if not isinstance(x, (list, tuple)): x = list(x) #raise TypeError, "x must be a list or tuple" self.__hash = None self.__cr = cr if cr_str is None: self.__cr_str = cr else: self.__cr_str = cr_str if isinstance(x, Sequence): if universe is None or universe == x.__universe: list.__init__(self, x) self.__universe = x.__universe self._is_immutable = immutable return if universe is None: if len(x) == 0: import sage.categories.all universe = sage.categories.all.Objects() else: import sage.structure.element as coerce y = x x = list( x ) # make a copy, or we'd change the type of the elements of x, which would be bad. if use_sage_types: # convert any Python builtin numerical types to Sage objects from sage.rings.integer_ring import ZZ from sage.rings.real_double import RDF from sage.rings.complex_double import CDF for i in range(len(x)): if isinstance(x[i], int) or isinstance(x[i], long): x[i] = ZZ(x[i]) elif isinstance(x[i], float): x[i] = RDF(x[i]) elif isinstance(x[i], complex): x[i] = CDF(x[i]) # start the pairwise coercion for i in range(len(x) - 1): try: x[i], x[i + 1] = coerce.canonical_coercion( x[i], x[i + 1]) except TypeError: import sage.categories.all universe = sage.categories.all.Objects() x = list(y) check = False # no point break if universe is None: # no type errors raised. universe = coerce.parent(x[len(x) - 1]) #universe = sage.structure.coerce.parent(x[0]) self.__universe = universe if check: x = [universe(t) for t in x] list.__init__(self, x) self._is_immutable = immutable
def numerical_approx(x, prec=None, digits=None, algorithm=None): r""" Returns a numerical approximation of an object ``x`` with at least ``prec`` bits (or decimal ``digits``) of precision. .. note:: Both upper case ``N`` and lower case ``n`` are aliases for :func:`numerical_approx`, and all three may be used as methods. INPUT: - ``x`` - an object that has a numerical_approx method, or can be coerced into a real or complex field - ``prec`` (optional) - an integer (bits of precision) - ``digits`` (optional) - an integer (digits of precision) - ``algorithm`` (optional) - a string specifying the algorithm to use for functions that implement more than one If neither the ``prec`` or ``digits`` are specified, the default is 53 bits of precision. If both are specified, then ``prec`` is used. EXAMPLES:: sage: numerical_approx(pi, 10) 3.1 sage: numerical_approx(pi, digits=10) 3.141592654 sage: numerical_approx(pi^2 + e, digits=20) 12.587886229548403854 sage: n(pi^2 + e) 12.5878862295484 sage: N(pi^2 + e) 12.5878862295484 sage: n(pi^2 + e, digits=50) 12.587886229548403854194778471228813633070946500941 sage: a = CC(-5).n(prec=40) sage: b = ComplexField(40)(-5) sage: a == b True sage: parent(a) is parent(b) True sage: numerical_approx(9) 9.00000000000000 You can also usually use method notation. :: sage: (pi^2 + e).n() 12.5878862295484 sage: (pi^2 + e).N() 12.5878862295484 sage: (pi^2 + e).numerical_approx() 12.5878862295484 Vectors and matrices may also have their entries approximated. :: sage: v = vector(RDF, [1,2,3]) sage: v.n() (1.00000000000000, 2.00000000000000, 3.00000000000000) sage: v = vector(CDF, [1,2,3]) sage: v.n() (1.00000000000000, 2.00000000000000, 3.00000000000000) sage: _.parent() Vector space of dimension 3 over Complex Field with 53 bits of precision sage: v.n(prec=20) (1.0000, 2.0000, 3.0000) sage: u = vector(QQ, [1/2, 1/3, 1/4]) sage: n(u, prec=15) (0.5000, 0.3333, 0.2500) sage: n(u, digits=5) (0.50000, 0.33333, 0.25000) sage: v = vector(QQ, [1/2, 0, 0, 1/3, 0, 0, 0, 1/4], sparse=True) sage: u = v.numerical_approx(digits=4) sage: u.is_sparse() True sage: u (0.5000, 0.0000, 0.0000, 0.3333, 0.0000, 0.0000, 0.0000, 0.2500) sage: A = matrix(QQ, 2, 3, range(6)) sage: A.n() [0.000000000000000 1.00000000000000 2.00000000000000] [ 3.00000000000000 4.00000000000000 5.00000000000000] sage: B = matrix(Integers(12), 3, 8, srange(24)) sage: N(B, digits=2) [0.00 1.0 2.0 3.0 4.0 5.0 6.0 7.0] [ 8.0 9.0 10. 11. 0.00 1.0 2.0 3.0] [ 4.0 5.0 6.0 7.0 8.0 9.0 10. 11.] Internally, numerical approximations of real numbers are stored in base-2. Therefore, numbers which look the same in their decimal expansion might be different:: sage: x=N(pi, digits=3); x 3.14 sage: y=N(3.14, digits=3); y 3.14 sage: x==y False sage: x.str(base=2) '11.001001000100' sage: y.str(base=2) '11.001000111101' Increasing the precision of a floating point number is not allowed:: sage: CC(-5).n(prec=100) Traceback (most recent call last): ... TypeError: cannot approximate to a precision of 100 bits, use at most 53 bits sage: n(1.3r, digits=20) Traceback (most recent call last): ... TypeError: cannot approximate to a precision of 70 bits, use at most 53 bits sage: RealField(24).pi().n() Traceback (most recent call last): ... TypeError: cannot approximate to a precision of 53 bits, use at most 24 bits As an exceptional case, ``digits=1`` usually leads to 2 digits (one significant) in the decimal output (see :trac:`11647`):: sage: N(pi, digits=1) 3.2 sage: N(pi, digits=2) 3.1 sage: N(100*pi, digits=1) 320. sage: N(100*pi, digits=2) 310. In the following example, ``pi`` and ``3`` are both approximated to two bits of precision and then subtracted, which kills two bits of precision:: sage: N(pi, prec=2) 3.0 sage: N(3, prec=2) 3.0 sage: N(pi - 3, prec=2) 0.00 TESTS:: sage: numerical_approx(I) 1.00000000000000*I sage: x = QQ['x'].gen() sage: F.<k> = NumberField(x^2+2, embedding=sqrt(CC(2))*CC.0) sage: numerical_approx(k) 1.41421356237309*I sage: type(numerical_approx(CC(1/2))) <type 'sage.rings.complex_number.ComplexNumber'> The following tests :trac:`10761`, in which ``n()`` would break when called on complex-valued algebraic numbers. :: sage: E = matrix(3, [3,1,6,5,2,9,7,3,13]).eigenvalues(); E [18.16815365088822?, -0.08407682544410650? - 0.2190261484802906?*I, -0.08407682544410650? + 0.2190261484802906?*I] sage: E[1].parent() Algebraic Field sage: [a.n() for a in E] [18.1681536508882, -0.0840768254441065 - 0.219026148480291*I, -0.0840768254441065 + 0.219026148480291*I] Make sure we've rounded up log(10,2) enough to guarantee sufficient precision (:trac:`10164`):: sage: ks = 4*10**5, 10**6 sage: check_str_length = lambda k: len(str(numerical_approx(1+10**-k,digits=k+1)))-1 >= k+1 sage: check_precision = lambda k: numerical_approx(1+10**-k,digits=k+1)-1 > 0 sage: all(check_str_length(k) and check_precision(k) for k in ks) True Testing we have sufficient precision for the golden ratio (:trac:`12163`), note that the decimal point adds 1 to the string length:: sage: len(str(n(golden_ratio, digits=5000))) 5001 sage: len(str(n(golden_ratio, digits=5000000))) # long time (4s on sage.math, 2012) 5000001 Check that :trac:`14778` is fixed:: sage: n(0, algorithm='foo') 0.000000000000000 """ if prec is None: if digits is None: prec = 53 else: prec = int((digits+1) * LOG_TEN_TWO_PLUS_EPSILON) + 1 try: return x._numerical_approx(prec, algorithm=algorithm) except AttributeError: pass from sage.structure.element import parent B = parent(x) RR = sage.rings.real_mpfr.RealField(prec) map = RR.coerce_map_from(B) if map is not None: return map(x) CC = sage.rings.complex_field.ComplexField(prec) map = CC.coerce_map_from(B) if map is not None: return map(x) # Coercion didn't work: there are 3 possibilities: # (1) There is a coercion possible to a lower precision # (2) There is a conversion but no coercion # (3) The type doesn't convert at all # Figure out input precision to check for case (1) try: inprec = x.prec() except AttributeError: if prec > 53 and CDF.has_coerce_map_from(B): # If we can coerce to CDF, assume input precision was 53 bits inprec = 53 else: # Otherwise, assume precision wasn't the issue inprec = prec if prec > inprec: raise TypeError("cannot approximate to a precision of %s bits, use at most %s bits" % (prec, inprec)) # The issue is not precision, try conversion instead try: return RR(x) except (TypeError, ValueError): pass return CC(x)
def Zinf(a, l, m, r, algorithm='spline'): r""" Amplitude factor of the mode `(\ell,m)`. The factor `Z^\infty_{\ell m}(r)` is obtained by spline interpolation of tabulated numerical solutions of the radial component of the Teukolsky equation. INPUT: - ``a`` -- BH angular momentum parameter (in units of `M`, the BH mass) - ``l`` -- integer >= 2; the harmonic degree `\ell` - ``m`` -- integer within the range ``[-l, l]``; the azimuthal number `m` - ``r`` -- areal radius of the orbit (in units of `M`) - ``algorithm`` -- (default: ``'spline'``) string describing the computational method; allowed values are - ``'spline'``: spline interpolation of tabulated data - ``'1.5PN'`` (only for ``a=0``): 1.5-post-Newtonian expansion following E. Poisson, Phys. Rev. D **47**, 1497 (1993) [:doi:`10.1103/PhysRevD.47.1497`], with a minus one factor accounting for a different convention for the metric signature. OUTPUT: - coefficient `Z^\infty_{\ell m}(r)` (in units of `M^{-2}`) EXAMPLES:: sage: from kerrgeodesic_gw import Zinf sage: Zinf(0.98, 2, 2, 1.7) # tol 1.0e-13 -0.04302234478778856 + 0.28535368610053824*I sage: Zinf(0., 2, 2, 10.) # tol 1.0e-13 0.0011206407919254163 - 0.0003057608384581628*I sage: Zinf(0., 2, 2, 10., algorithm='1.5PN') # tol 1.0e-13 0.0011971529546749354 - 0.0003551610880408921*I """ if m < 0: return (-1)**l * Zinf(a, l, -m, r, algorithm=algorithm).conjugate() if algorithm == '1.5PN': if a == 0: # the factor (-1) below accounts for a difference of signature with # Poisson (1993): return -Zinf_Schwarzchild_PN(l, m, r) raise ValueError("a must be zero for algorithm='1.5PN'") a = RDF(a) param = (a, l, m) if param in _cached_splines: splines = _cached_splines[param] else: file_name = "data/Zinf_a{:.1f}.dat".format(float(a)) if a <= 0.9 \ else "data/Zinf_a{:.2f}.dat".format(float(a)) file_name = os.path.join(os.path.dirname(__file__), file_name) r_high_l = 20. if a <= 0.9 else 10. with open(file_name, "r") as data_file: lm_values = [] # l values up to 10 (for r <= r_high_l) lm_values_low = [] # l values up to 5 only (for r > r_high_l) for ld in range(2, 6): for md in range(1, ld + 1): lm_values_low.append((ld, md)) for ld in range(2, 11): for md in range(1, ld + 1): lm_values.append((ld, md)) Zreal = {} Zimag = {} for (ld, md) in lm_values: Zreal[(ld, md)] = [] Zimag[(ld, md)] = [] for line in data_file: items = line.split('\t') rd = RDF(items.pop(0)) if rd <= r_high_l: for (ld, md) in lm_values: Zreal[(ld, md)].append((rd, RDF(items.pop(0)))) Zimag[(ld, md)].append((rd, RDF(items.pop(0)))) else: for (ld, md) in lm_values_low: Zreal[(ld, md)].append((rd, RDF(items.pop(0)))) Zimag[(ld, md)].append((rd, RDF(items.pop(0)))) for (ld, md) in lm_values: _cached_splines[(a, ) + (ld, md)] = (spline(Zreal[(ld, md)]), spline(Zimag[(ld, md)])) if param not in _cached_splines: raise ValueError( "Zinf: case (a, l, m) = {} not implemented".format(param)) splines = _cached_splines[param] # The factor (-1)**(l+m) below accounts for a difference of convention # in the C++ code used to produce the data files return (-1)**(l + m) * CDF(splines[0](r), splines[1](r))
def Zinf_Schwarzchild_PN(l, m, r): r""" Amplitude factor of the mode `(\ell,m)` for a Schwarzschild BH at the 1.5PN level. The 1.5PN formulas are taken from E. Poisson, Phys. Rev. D **47**, 1497 (1993), :doi:`10.1103/PhysRevD.47.1497`. INPUT: - ``l`` -- integer >= 2; the harmonic degree `\ell` - ``m`` -- integer within the range ``[-l, l]``; the azimuthal number `m` - ``r`` -- areal radius of the orbit (in units of `M`, the BH mass) OUTPUT: - coefficient `Z^\infty_{\ell m}(r)` (in units of `M^{-2}`) EXAMPLES:: sage: from kerrgeodesic_gw import Zinf_Schwarzchild_PN sage: Zinf_Schwarzchild_PN(2, 2, 6.) # tol 1.0e-13 -0.00981450418730346 + 0.003855681972781947*I sage: Zinf_Schwarzchild_PN(5, 3, 6.) # tol 1.0e-13 -6.958527913913504e-05*I """ if m < 0: return (-1)**l * Zinf_Schwarzchild_PN(l, -m, r).conjugate() v = 1. / sqrt(r) if l == 2: b = RDF(sqrt(pi / 5.) / r**4) if m == 1: return CDF(4. / 3. * I * b * v * (1 - 17. / 28. * v**2)) if m == 2: return CDF(-16 * b * (1 - 107. / 42. * v**2 + (2 * pi + 4 * I * (3 * ln(2 * v) - 0.839451001765134)) * v**3)) if l == 3: b = RDF(sqrt(pi / 7.) / r**(4.5)) if m == 1: return CDF(I / 3. * b / sqrt(10.) * (1 - 8. / 3. * v**2)) if m == 2: return CDF(-16. / 3. * b * v) if m == 3: return CDF(-81 * I * b / sqrt(8.) * (1 - 4 * v**2)) if l == 4: b = RDF(sqrt(pi) / r**5) if m == 1: return CDF(I / 105. * b / sqrt(2) * v) if m == 2: return CDF(-16. / 63. * b) if m == 3: return CDF(-81. / 5. * I * b / sqrt(14.) * v) if m == 4: return CDF(512. / 9. * b / sqrt(7.)) if l == 5: b = RDF(sqrt(pi) / r**(5.5)) if m == 1: return CDF(I / 360. * b / sqrt(77.)) if m == 2: return CDF(0) if m == 3: return CDF(-81. / 40. * I * b * sqrt(3. / 22.)) if m == 4: return CDF(0) if m == 5: return CDF(3125. / 24. * I * b * sqrt(5. / 66.)) raise NotImplementedError("{} not implemented".format((l, m)))
def Sequence(x, universe=None, check=True, immutable=False, cr=False, cr_str=None, use_sage_types=False): """ A mutable list of elements with a common guaranteed universe, which can be set immutable. A universe is either an object that supports coercion (e.g., a parent), or a category. INPUT: - ``x`` - a list or tuple instance - ``universe`` - (default: None) the universe of elements; if None determined using canonical coercions and the entire list of elements. If list is empty, is category Objects() of all objects. - ``check`` -- (default: True) whether to coerce the elements of x into the universe - ``immutable`` - (default: True) whether or not this sequence is immutable - ``cr`` - (default: False) if True, then print a carriage return after each comma when printing this sequence. - ``cr_str`` - (default: False) if True, then print a carriage return after each comma when calling ``str()`` on this sequence. - ``use_sage_types`` -- (default: False) if True, coerce the built-in Python numerical types int, long, float, complex to the corresponding Sage types (this makes functions like vector() more flexible) OUTPUT: - a sequence EXAMPLES:: sage: v = Sequence(range(10)) sage: v.universe() <type 'int'> sage: v [0, 1, 2, 3, 4, 5, 6, 7, 8, 9] We can request that the built-in Python numerical types be coerced to Sage objects:: sage: v = Sequence(range(10), use_sage_types=True) sage: v.universe() Integer Ring sage: v [0, 1, 2, 3, 4, 5, 6, 7, 8, 9] You can also use seq for "Sequence", which is identical to using Sequence:: sage: v = seq([1,2,1/1]); v [1, 2, 1] sage: v.universe() Rational Field sage: v.parent() Category of sequences in Rational Field sage: v.parent()([3,4/3]) [3, 4/3] Note that assignment coerces if possible,:: sage: v = Sequence(range(10), ZZ) sage: a = QQ(5) sage: v[3] = a sage: parent(v[3]) Integer Ring sage: parent(a) Rational Field sage: v[3] = 2/3 Traceback (most recent call last): ... TypeError: no conversion of this rational to integer Sequences can be used absolutely anywhere lists or tuples can be used:: sage: isinstance(v, list) True Sequence can be immutable, so entries can't be changed:: sage: v = Sequence([1,2,3], immutable=True) sage: v.is_immutable() True sage: v[0] = 5 Traceback (most recent call last): ... ValueError: object is immutable; please change a copy instead. Only immutable sequences are hashable (unlike Python lists), though the hashing is potentially slow, since it first involves conversion of the sequence to a tuple, and returning the hash of that.:: sage: v = Sequence(range(10), ZZ, immutable=True) sage: hash(v) 1591723448 # 32-bit -4181190870548101704 # 64-bit If you really know what you are doing, you can circumvent the type checking (for an efficiency gain):: sage: list.__setitem__(v, int(1), 2/3) # bad circumvention sage: v [0, 2/3, 2, 3, 4, 5, 6, 7, 8, 9] sage: list.__setitem__(v, int(1), int(2)) # not so bad circumvention You can make a sequence with a new universe from an old sequence.:: sage: w = Sequence(v, QQ) sage: w [0, 2, 2, 3, 4, 5, 6, 7, 8, 9] sage: w.universe() Rational Field sage: w[1] = 2/3 sage: w [0, 2/3, 2, 3, 4, 5, 6, 7, 8, 9] Sequences themselves live in a category, the category of all sequences in the given universe.:: sage: w.category() Category of sequences in Rational Field This is also the parent of any sequence:: sage: w.parent() Category of sequences in Rational Field The default universe for any sequence, if no compatible parent structure can be found, is the universe of all Sage objects. This example illustrates how every element of a list is taken into account when constructing a sequence.:: sage: v = Sequence([1,7,6,GF(5)(3)]); v [1, 2, 1, 3] sage: v.universe() Finite Field of size 5 sage: v.parent() Category of sequences in Finite Field of size 5 sage: v.parent()([7,8,9]) [2, 3, 4] """ from sage.rings.polynomial.multi_polynomial_ideal import MPolynomialIdeal if isinstance(x, Sequence_generic) and universe is None: universe = x.universe() x = list(x) if isinstance(x, MPolynomialIdeal) and universe is None: universe = x.ring() x = x.gens() if universe is None: if not isinstance(x, (list, tuple)): x = list(x) #raise TypeError("x must be a list or tuple") if len(x) == 0: import sage.categories.all universe = sage.categories.all.Objects() else: import sage.structure.element as coerce y = x x = list(x) # make a copy, or we'd change the type of the elements of x, which would be bad. if use_sage_types: # convert any Python built-in numerical types to Sage objects from sage.rings.integer_ring import ZZ from sage.rings.real_double import RDF from sage.rings.complex_double import CDF for i in range(len(x)): if isinstance(x[i], int) or isinstance(x[i], long): x[i] = ZZ(x[i]) elif isinstance(x[i], float): x[i] = RDF(x[i]) elif isinstance(x[i], complex): x[i] = CDF(x[i]) # start the pairwise coercion for i in range(len(x)-1): try: x[i], x[i+1] = coerce.canonical_coercion(x[i],x[i+1]) except TypeError: import sage.categories.all universe = sage.categories.all.Objects() x = list(y) check = False # no point break if universe is None: # no type errors raised. universe = coerce.parent(x[len(x)-1]) from sage.rings.polynomial.multi_polynomial_ring import is_MPolynomialRing from sage.rings.quotient_ring import is_QuotientRing from sage.rings.polynomial.pbori import BooleanMonomialMonoid if is_MPolynomialRing(universe) or \ (is_QuotientRing(universe) and is_MPolynomialRing(universe.cover_ring())) or \ isinstance(universe, BooleanMonomialMonoid): from sage.rings.polynomial.multi_polynomial_sequence import PolynomialSequence try: return PolynomialSequence(x, universe, immutable=immutable, cr=cr, cr_str=cr_str) except (TypeError,AttributeError): return Sequence_generic(x, universe, check, immutable, cr, cr_str, use_sage_types) else: return Sequence_generic(x, universe, check, immutable, cr, cr_str, use_sage_types)
def julia_plot(f=None, **kwds): r""" Plots the Julia set of a given polynomial ``f``. Users can specify whether they would like to display the Mandelbrot side by side with the Julia set with the ``mandelbrot`` argument. If ``f`` is not specified, this method defaults to `f(z) = z^2-1`. The Julia set of a polynomial ``f`` is the set of complex numbers `z` for which the function `f(z)` is bounded under iteration. The Julia set can be visualized by plotting each point in the set in the complex plane. Julia sets are examples of fractals when plotted in the complex plane. ALGORITHM: Let `R_c = \bigl(1 + \sqrt{1 + 4|c|}\bigr)/2` if the polynomial is of the form `f(z) = z^2 + c`; otherwise, let `R_c = 2`. For every `p \in \mathbb{C}`, if `|f^{k}(p)| > R_c` for some `k \geq 0`, then `f^{n}(p) \to \infty`. Let `N` be the maximum number of iterations. Compute the first `N` points on the orbit of `p` under `f`. If for any `k < N`, `|f^{k}(p)| > R_c`, we stop the iteration and assign a color to the point `p` based on how quickly `p` escaped to infinity under iteration of `f`. If `|f^{i}(p)| \leq R_c` for all `i \leq N`, we assume `p` is in the Julia set and assign the point `p` the color black. INPUT: - ``f`` -- input polynomial (optional - default: ``z^2 - 1``). - ``period`` -- list (optional - default: ``None``), returns the Julia set for a random `c` value with the given (formal) cycle structure. - ``mandelbrot`` -- boolean (optional - default: ``True``), when set to ``True``, an image of the Mandelbrot set is appended to the right of the Julia set. - ``point_color`` -- RGB color (optional - default: ``'tomato'``), color of the point `c` in the Mandelbrot set (any valid input for Color). - ``x_center`` -- double (optional - default: ``-1.0``), Real part of center point. - ``y_center`` -- double (optional - default: ``0.0``), Imaginary part of center point. - ``image_width`` -- double (optional - default: ``4.0``), width of image in the complex plane. - ``max_iteration`` -- long (optional - default: ``500``), maximum number of iterations the map `f(z)`. - ``pixel_count`` -- long (optional - default: ``500``), side length of image in number of pixels. - ``base_color`` -- hex color (optional - default: ``'steelblue'``), color used to determine the coloring of set (any valid input for Color). - ``level_sep`` -- long (optional - default: 1), number of iterations between each color level. - ``number_of_colors`` -- long (optional - default: 30), number of colors used to plot image. - ``interact`` -- boolean (optional - default: ``False``), controls whether plot will have interactive functionality. OUTPUT: 24-bit RGB image of the Julia set in the complex plane. .. TODO:: Implement the side-by-side Mandelbrot-Julia plots for general one-parameter families of polynomials. EXAMPLES: The default ``f`` is `z^2 - 1`:: sage: julia_plot() 1001x500px 24-bit RGB image To display only the Julia set, set ``mandelbrot`` to ``False``:: sage: julia_plot(mandelbrot=False) 500x500px 24-bit RGB image :: sage: R.<z> = CC[] sage: f = z^3 - z + 1 sage: julia_plot(f) 500x500px 24-bit RGB image To display an interactive plot of the Julia set in the Notebook, set ``interact`` to ``True``. (This is only implemented for polynomials of the form ``f = z^2 + c``):: sage: julia_plot(interact=True) interactive(children=(FloatSlider(value=-1.0, description='Real c'... :: sage: R.<z> = CC[] sage: f = z^2 + 1/2 sage: julia_plot(f,interact=True) interactive(children=(FloatSlider(value=0.5, description='Real c'... To return the Julia set of a random `c` value with (formal) cycle structure `(2,3)`, set ``period = [2,3]``:: sage: julia_plot(period=[2,3]) 1001x500px 24-bit RGB image To return all of the Julia sets of `c` values with (formal) cycle structure `(2,3)`:: sage: period = [2,3] # not tested ....: R.<c> = QQ[] ....: P.<x,y> = ProjectiveSpace(R,1) ....: f = DynamicalSystem([x^2+c*y^2, y^2]) ....: L = f.dynatomic_polynomial(period).subs({x:0,y:1}).roots(ring=CC) ....: c_values = [k[0] for k in L] ....: for c in c_values: ....: julia_plot(c) Polynomial maps can be defined over a polynomial ring or a fraction field, so long as ``f`` is polynomial:: sage: R.<z> = CC[] sage: f = z^2 - 1 sage: julia_plot(f) 1001x500px 24-bit RGB image :: sage: R.<z> = CC[] sage: K = R.fraction_field(); z = K.gen() sage: f = z^2-1 sage: julia_plot(f) 1001x500px 24-bit RGB image Interact functionality is not implemented if the polynomial is not of the form `f = z^2 + c`:: sage: R.<z> = CC[] sage: f = z^3 + 1 sage: julia_plot(f, interact=True) Traceback (most recent call last): ... NotImplementedError: The interactive plot is only implemented for ... """ # extract keyword arguments period = kwds.pop("period", None) mandelbrot = kwds.pop("mandelbrot", True) point_color = kwds.pop("point_color", 'tomato') x_center = kwds.pop("x_center", 0.0) y_center = kwds.pop("y_center", 0.0) image_width = kwds.pop("image_width", 4.0) max_iteration = kwds.pop("max_iteration", 500) pixel_count = kwds.pop("pixel_count", 500) base_color = kwds.pop("base_color", 'steelblue') level_sep = kwds.pop("level_sep", 1) number_of_colors = kwds.pop("number_of_colors", 30) interacts = kwds.pop("interact", False) f_is_default_after_all = None if period: # pick a random c with the specified period R = PolynomialRing(CC, 'c') c = R.gen() x, y = ProjectiveSpace(R, 1, 'x,y').gens() F = DynamicalSystem([x**2 + c * y**2, y**2]) L = F.dynatomic_polynomial(period).subs({x: 0, y: 1}).roots(ring=CC) c = L[randint(0, len(L) - 1)][0] base_color = Color(base_color) point_color = Color(point_color) EPS = 0.00001 if f is not None and period is None: # f user-specified and no period given # try to coerce f to live in a polynomial ring S = PolynomialRing(CC, names='z') z = S.gen() try: f_poly = S(f) except TypeError: R = f.parent() if not (R.is_integral_domain() and (CC.is_subring(R) or CDF.is_subring(R))): raise ValueError('Given `f` must be a complex polynomial.') else: raise NotImplementedError( 'Julia sets not implemented for rational functions.' ) if (f_poly - z*z) in CC: # f is specified and of the form z^2 + c. f_is_default_after_all = True c = f_poly - z*z else: # f is specified and not of the form z^2 + c if interacts: raise NotImplementedError( "The interactive plot is only implemented for " "polynomials of the form f = z^2 + c." ) else: return general_julia(f_poly, x_center, y_center, image_width, max_iteration, pixel_count, level_sep, number_of_colors, base_color) # otherwise we can use fast_julia_plot for z^2 + c if f_is_default_after_all or f is None or period is not None: # specify default c = -1 value if f and period were not specified if not f_is_default_after_all and period is None: c = -1 c = CC(c) c_real = c.real() c_imag = c.imag() if interacts: # set widgets from ipywidgets.widgets import FloatSlider, IntSlider, \ ColorPicker, interact widgets = dict( c_real = FloatSlider(min=-2.0, max=2.0, step=EPS, value=c_real, description="Real c"), c_imag = FloatSlider(min=-2.0, max=2.0, step=EPS, value=c_imag, description="Imag c"), x_center = FloatSlider(min=-1.0, max=1.0, step=EPS, value=x_center, description="Real center"), y_center = FloatSlider(min=-1.0, max=1.0, step=EPS, value=y_center, description="Imag center"), image_width = FloatSlider(min=EPS, max=4.0, step=EPS, value=image_width, description="Width"), max_iteration = IntSlider(min=0, max=1000, value=max_iteration, description="Iterations"), pixel_count = IntSlider(min=10, max=1000, value=pixel_count, description="Pixels"), level_sep = IntSlider(min=1, max=20, value=level_sep, description="Color sep"), color_num = IntSlider(min=1, max=100, value=number_of_colors, description="# Colors"), base_color = ColorPicker(value=base_color.html_color(), description="Base color"), ) if mandelbrot: widgets["point_color"] = ColorPicker(value=point_color.html_color(), description="Point color") return interact(**widgets).widget(julia_helper) else: return interact(**widgets).widget(fast_julia_plot) elif mandelbrot: # non-interactive with mandelbrot return julia_helper(c_real, c_imag, x_center, y_center, image_width, max_iteration, pixel_count, level_sep, number_of_colors, base_color, point_color) else: # non-interactive without mandelbrot return fast_julia_plot(c_real, c_imag, x_center, y_center, image_width, max_iteration, pixel_count, level_sep, number_of_colors, base_color)
def h_toy_model_semi_analytic(u, theta, phi, a, r0, phi0, lam, Dphi, l_max=10): r""" Return the gravitational wave emitted by a matter blob orbiting a Kerr black hole (semi-analytic computation based on a toy model surface density). The surface density of the matter blob is that given by :func:`surface_density_toy_model`. The gravitational wave is computed according to the formula .. MATH:: h = \frac{2\mu}{r} \, \sum_{\ell=2}^{\infty} \sum_{m=-\ell}^\ell \frac{Z^\infty_{\ell m}(r_0)}{(m\omega_0)^2} \; \text{sinc}\left( \frac{m}{2} \Delta\varphi \right) \, \text{sinc}\left( \frac{3}{4} \varepsilon \, m \omega_0 (1-a\omega_0)u \right) e^{- i m (\omega_0 u + \phi_0)} \, _{-2}S_{\ell m}^{a m \omega_0}(\theta,\varphi) INPUT: - ``u`` -- retarded time coordinate of the observer (in units of `M`, the BH mass): `u = t - r_*`, where `t` is the Boyer-Lindquist time coordinate and `r_*` is the tortoise coordinate - ``theta`` -- Boyer-Lindquist colatitute `\theta` of the observer - ``phi`` -- Boyer-Lindquist azimuthal coordinate `\phi` of the observer - ``a`` -- BH angular momentum parameter (in units of `M`) - ``r0`` -- mean radius `r_0` of the matter blob (Boyer-Lindquist coordinate) - ``phi0`` -- mean azimuthal angle `\phi_0` of the matter blob (Boyer-Lindquist coordinate) - ``lam`` -- radial extent `\lambda` of the matter blob - ``Dphi``-- opening angle `\Delta\phi` of the matter blob - ``l_max`` -- (default: 10) upper bound in the summation over the harmonic degree `\ell` OUTPUT: - a pair ``(hp, hc)``, where ``hp`` (resp. ``hc``) is `(r / \mu) h_+` (resp. `(r / \mu) h_\times`), `\mu` being the blob's mass and `r` is the Boyer-Lindquist radial coordinate of the observer EXAMPLES: Schwarzschild black hole:: sage: from kerrgeodesic_gw import h_toy_model_semi_analytic sage: a = 0 sage: r0, phi0, lam, Dphi = 6.5, 0, 0.6, 0.1 sage: u = 60. sage: h_toy_model_semi_analytic(u, pi/4, 0., a, r0, phi0, lam, Dphi) # tol 1.0e-13 (0.2999183296797872, 0.36916647790743246) sage: hp, hc = _ Comparison with the exact value:: sage: from kerrgeodesic_gw import (h_blob, blob_mass, ....: surface_density_toy_model) sage: param_surf_dens = [r0, phi0, lam, Dphi] sage: integ_range = [6.2, 6.8, -0.05, 0.05] sage: mu = blob_mass(a, surface_density_toy_model, param_surf_dens, ....: integ_range)[0] sage: hp0 = h_blob(u, pi/4, 0., a, surface_density_toy_model, ....: param_surf_dens, integ_range)[0] / mu sage: hc0 = h_blob(u, pi/4, 0., a, surface_density_toy_model, ....: param_surf_dens, integ_range, mode='x')[0] / mu sage: hp0, hc0 # tol 1.0e-13 (0.2951163078053617, 0.3743683023327848) sage: (hp - hp0) / hp0 # tol 1.0e-13 0.01627162494047128 sage: (hc - hc0) / hc0 # tol 1.0e-13 -0.013894938201066784 """ import numpy from sage.rings.real_double import RDF from sage.rings.complex_double import CDF from sage.symbolic.all import i as I from .spin_weighted_spherical_harm import spin_weighted_spherical_harmonic from .spin_weighted_spheroidal_harm import spin_weighted_spheroidal_harmonic from .zinf import Zinf u = RDF(u) theta = RDF(theta) phi = RDF(phi) a = RDF(a) omega0 = RDF(1. / (r0**1.5 + a)) eps = lam/r0 resu = CDF(0) for l in range(2, l_max+1): for m in range(-l, l+1): if m == 0: # m=0 is skipped continue # m_omega0 = RDF(m*omega0) if a == 0: Slm = spin_weighted_spherical_harmonic(-2, l, m, theta, phi, numerical=RDF) else: a = RDF(a) Slm = spin_weighted_spheroidal_harmonic(-2, l, m, a*m_omega0, theta, phi) # Division by pi in the Sinc function due to the defintion used by numpy resu += Zinf(a, l, m, r0) / m_omega0**2 \ * numpy.sinc(m*Dphi/2./numpy.pi) \ * numpy.sinc(0.75*eps*m_omega0*(1-a*omega0)*u/numpy.pi) \ * CDF(exp(-I*(m_omega0*u + m*phi0))) * Slm resu *= 2 return (resu.real(), -resu.imag())