def approx_abs_real(self, prec): r""" Compute an approximation with absolute error about 2^(-prec). """ if isinstance(self.value.parent(), RealBallField): return self.value elif self.value.is_zero(): return RealBallField(max(2, prec)).zero() elif self.is_real(): expo = ZZ(IR(self.value).abs().log(2).upper().ceil()) rel_prec = max(2, prec + expo + 10) val = RealBallField(rel_prec)(self.value) return val else: raise ValueError("point may not be real")
def _test_fun_approx(pol, ref, disk_rad=None, interval_rad=None, prec=53, test_count=100): r""" EXAMPLES:: sage: from ore_algebra.analytic.polynomial_approximation import _test_fun_approx sage: _test_fun_approx(lambda x: x.exp(), lambda x: x.exp() + x/1000, ....: interval_rad=1) Traceback (most recent call last): ... AssertionError: z = ..., ref(z) = ... not in pol(z) = ... """ from sage.rings.real_mpfr import RealField from sage.rings.real_arb import RealBallField from sage.rings.complex_arb import ComplexBallField my_RR = RealField(prec) my_RBF = RealBallField(prec) my_CBF = ComplexBallField(prec) if bool(disk_rad) == bool(interval_rad): raise ValueError rad = disk_rad or interval_rad for _ in range(test_count): rho = my_RBF(my_RR.random_element(-rad, rad)) if disk_rad: exp_i_theta = my_CBF(my_RR.random_element(0, 1)).exppii() z = rho*exp_i_theta elif interval_rad: z = rho ref_z = ref(z) pol_z = pol(z) if not ref_z.overlaps(pol_z): fmt = "z = {}, ref(z) = {} not in pol(z) = {}" raise AssertionError(fmt.format(z, ref_z, pol_z))
def analytic_continuation(ctx, ini=None, post=None): """ INPUT: - ``ini`` (constant matrix, optional) - initial values, one column per solution - ``post`` (matrix of polynomial/rational functions, optional) - linear combinations of the first Taylor coefficients to take, as a function of the evaluation point TESTS:: sage: from ore_algebra import DifferentialOperators sage: _, x, Dx = DifferentialOperators() sage: (Dx^2 + 2*x*Dx).numerical_solution([0, 2/sqrt(pi)], [0,i]) [+/- ...] + [1.65042575879754...]*I """ logger.info("path: %s", ctx.path) eps1 = (ctx.eps / (1 + len(ctx.path))) >> 2 # TBI, +: move to ctx? prec = utilities.prec_from_eps(eps1) if ini is not None: if not isinstance(ini, Matrix): # should this be here? try: ini = matrix(ctx.dop.order(), 1, list(ini)) except (TypeError, ValueError): raise ValueError("incorrect initial values: {}".format(ini)) try: ini = ini.change_ring(RealBallField(prec)) except (TypeError, ValueError): ini = ini.change_ring(ComplexBallField(prec)) res = [] path_mat = identity_matrix(ZZ, ctx.dop.order()) def store_value_if_wanted(point): if point.options.get('keep_value'): value = path_mat if ini is not None: value = value * ini if post is not None: value = post(point.value) * value res.append((point.value, value)) store_value_if_wanted(ctx.path.vert[0]) for step in ctx.path: step_mat = step_transition_matrix(step, eps1, ctx=ctx) path_mat = step_mat * path_mat store_value_if_wanted(step.end) cm = sage.structure.element.get_coercion_model() OutputIntervals = cm.common_parent( utilities.ball_field(ctx.eps, ctx.real()), *[mat.base_ring() for pt, mat in res]) return [(pt, mat.change_ring(OutputIntervals)) for pt, mat in res]
def _disk(self, pt): assert pt.is_real() # Since approximation disks satisfy 2·rad ≤ dist(center, sing), any # approximation disk containing pt must have rad ≤ dist(pt, sing) max_rad = pt.dist_to_sing().min(self.max_rad) # What we want is the largest such disk containing pt expo = ZZ(max_rad.log(2).upper().ceil()) - 1 # rad = 2^expo logger.log(logging.DEBUG - 2, "max_rad = %s, expo = %s", max_rad, expo) while True: approx_pt = pt.approx_abs_real(-expo) mantissa = (approx_pt.squash() >> expo).floor() if ZZ(mantissa) % 2 == 0: mantissa += 1 center = mantissa << expo dist = Point(center, pt.dop).dist_to_sing() rad = RBF.one() << expo logger.log( logging.DEBUG - 2, "candidate disk: approx_pt = %s, mantissa = %s, " "center = %s, dist = %s, rad = %s", approx_pt, mantissa, center, dist, rad) if safe_ge(dist >> 1, rad): break expo -= 1 logger.debug("disk for %s: center=%s, rad=%s", pt, center, rad) # pt may be a ball with nonzero radius: check that it is contained in # our candidate disk log = RBF.zero() if 0 in approx_pt else approx_pt.abs().log(2) F = RealBallField(ZZ((expo - log).max(0).upper().ceil()) + 10) dist_to_center = (F(approx_pt) - F(center)).abs() if not safe_le(dist_to_center, rad): assert not safe_gt((approx_pt.squash() - center).squash(), rad) logger.info("check that |%s - %s| < %s failed", approx_pt, center, rad) return None, None # exactify center so that subsequent computations are not limited by the # precision of its parent center = QQ(center) return center, rad
def __init__(self, point, dop=None): """ TESTS:: sage: from ore_algebra import * sage: from ore_algebra.analytic.path import Point sage: Dops, x, Dx = DifferentialOperators() sage: [Point(z, Dx) ....: for z in [1, 1/2, 1+I, QQbar(I), RIF(1/3), CIF(1/3), pi, ....: RDF(1), CDF(I), 0.5r, 0.5jr, 10r, QQbar(1), AA(1/3)]] [1, 1/2, I + 1, I, [0.333333333333333...], [0.333333333333333...], 3.141592653589794?, 1.000000000000000, 1.000000000000000*I, 0.5000000000000000, 0.5000000000000000*I, 10, 1, 1/3] sage: Point(sqrt(2), Dx).iv() [1.414...] """ SageObject.__init__(self) from sage.rings.complex_double import ComplexDoubleField_class from sage.rings.complex_field import ComplexField_class from sage.rings.complex_interval_field import ComplexIntervalField_class from sage.rings.real_double import RealDoubleField_class from sage.rings.real_mpfi import RealIntervalField_class from sage.rings.real_mpfr import RealField_class point = sage.structure.coerce.py_scalar_to_element(point) try: parent = point.parent() except AttributeError: raise TypeError("unexpected value for point: " + repr(point)) if isinstance(point, Point): self.value = point.value elif isinstance( parent, (number_field_base.NumberField, RealBallField, ComplexBallField)): self.value = point elif QQ.has_coerce_map_from(parent): self.value = QQ.coerce(point) # must come before QQbar, due to a bogus coerce map (#14485) elif parent is sage.symbolic.ring.SR: try: return self.__init__(point.pyobject(), dop) except TypeError: pass try: return self.__init__(QQbar(point), dop) except (TypeError, ValueError, NotImplementedError): pass try: self.value = RLF(point) except (TypeError, ValueError): self.value = CLF(point) elif QQbar.has_coerce_map_from(parent): alg = QQbar.coerce(point) NF, val, hom = alg.as_number_field_element() if NF is QQ: self.value = QQ.coerce(val) # parent may be ZZ else: embNF = number_field.NumberField(NF.polynomial(), NF.variable_name(), embedding=hom(NF.gen())) self.value = val.polynomial()(embNF.gen()) elif isinstance( parent, (RealField_class, RealDoubleField_class, RealIntervalField_class)): self.value = RealBallField(point.prec())(point) elif isinstance(parent, (ComplexField_class, ComplexDoubleField_class, ComplexIntervalField_class)): self.value = ComplexBallField(point.prec())(point) else: try: self.value = RLF.coerce(point) except TypeError: self.value = CLF.coerce(point) parent = self.value.parent() assert (isinstance( parent, (number_field_base.NumberField, RealBallField, ComplexBallField)) or parent is RLF or parent is CLF) self.dop = dop or point.dop self.keep_value = False
def analytic_continuation(dop, path, eps, ctx=dctx, ini=None, post=None, return_local_bases=False): """ INPUT: - ``ini`` (constant matrix, optional) - initial values, one column per solution - ``post`` (matrix of polynomial/rational functions, optional) - linear combinations of the first Taylor coefficients to take, as a function of the evaluation point - ``return_local_bases`` (boolean) - if True, also compute and return the structure of local bases at all points where we are computing values of the solution OUTPUT: A list of dictionaries with information on the computed solution(s) at each evaluation point. TESTS:: sage: from ore_algebra import DifferentialOperators sage: _, x, Dx = DifferentialOperators() sage: (Dx^2 + 2*x*Dx).numerical_solution([0, 2/sqrt(pi)], [0,i]) [+/- ...] + [1.65042575879754...]*I """ if dop.is_zero(): raise ValueError("operator must be nonzero") _, _, _, dop = dop._normalize_base_ring() path = _process_path(dop, path, ctx) logger.info("path: %s", path) eps = bounds.IR(eps) eps1 = (eps/(1 + len(path))) >> 2 prec = utilities.prec_from_eps(eps1) if ini is not None: if not isinstance(ini, Matrix): # should this be here? try: ini = matrix(dop.order(), 1, list(ini)) except (TypeError, ValueError): raise ValueError("incorrect initial values: {}".format(ini)) try: ini = ini.change_ring(RealBallField(prec)) except (TypeError, ValueError): ini = ini.change_ring(ComplexBallField(prec)) def point_dict(point, value): if ini is not None: value = value*ini if post is not None and not post.is_one(): value = post(point.value)*value rec = {"point": point.value, "value": value} if return_local_bases: rec["structure"] = point.local_basis_structure() return rec res = [] z0 = path.vert[0] # XXX still imperfect in the case of a high-precision starting point with # relatively large radius... (do we care?) main = Step(z0, z0.simple_approx(ctx=ctx)) path_mat = step_transition_matrix(dop, main, eps1, ctx=ctx) if z0.keep_value(): res.append(point_dict(z0, identity_matrix(ZZ, dop.order()))) for step in path: main, dev = step.chain_simple(main.end, ctx=ctx) main_mat = step_transition_matrix(dop, main, eps1, ctx=ctx) path_mat = main_mat*path_mat if dev is not None: dev_mat = path_mat for sub in dev: sub_mat = step_transition_matrix(dop, sub, eps1, ctx=ctx) dev_mat = sub_mat*dev_mat res.append(point_dict(step.end, dev_mat)) cm = sage.structure.element.get_coercion_model() real = (rings.RIF.has_coerce_map_from(dop.base_ring().base_ring()) and all(v.is_real() for v in path.vert)) OutputIntervals = cm.common_parent( utilities.ball_field(eps, real), *[rec["value"].base_ring() for rec in res]) for rec in res: rec["value"] = rec["value"].change_ring(OutputIntervals) return res
def approx(self, pt, prec=None, post_transform=None): r""" TESTS:: sage: from ore_algebra import * sage: from ore_algebra.analytic.function import DFiniteFunction sage: DiffOps, x, Dx = DifferentialOperators() sage: h = DFiniteFunction(Dx^3-1, [0, 0, 1]) sage: h.approx(0, post_transform=Dx^2) [2.0000000000000...] sage: f = DFiniteFunction((x^2 + 1)*Dx^2 + 2*x*Dx, [0, 1], max_prec=20) sage: f.approx(1/3, prec=10) [0.32...] sage: f.approx(1/3, prec=40) [0.321750554396...] sage: f.approx(1/3, prec=10, post_transform=Dx) [0.9...] sage: f.approx(1/3, prec=40, post_transform=Dx) [0.900000000000...] sage: f.approx(1/3, prec=10, post_transform=Dx^2) [-0.54...] sage: f.approx(1/3, prec=40, post_transform=Dx^2) [-0.540000000000...] """ pt = Point(pt, self.dop) if prec is None: prec = _guess_prec(pt) if post_transform is None: post_transform = self.dop.parent().one() derivatives = min(post_transform.order() + 1, self._max_derivatives) post_transform = normalize_post_transform(self.dop, post_transform) if prec >= self.max_prec or not pt.is_real(): logger.info( "performing high-prec evaluation " "(pt=%s, prec=%s, post_transform=%s)", pt, prec, post_transform) ini, path = self._path_to(pt) eps = RBF.one() >> prec return self.dop.numerical_solution(ini, path, eps, post_transform=post_transform) center, rad = self._disk(pt) if center is None: # raise NotImplementedError logger.info("falling back on generic evaluator") ini, path = self._path_to(pt) eps = RBF.one() >> prec return self.dop.numerical_solution(ini, path, eps, post_transform=post_transform) approx = self._polys.get(center, []) Balls = RealBallField(prec) # due to the way the polynomials are recomputed, the precisions attached # to the successive derivatives are nonincreasing if (len(approx) < derivatives or approx[derivatives - 1].prec < prec): polys = self._update_approx(center, rad, prec, derivatives) else: polys = [a.pol for a in approx] bpt = Balls(pt.value) reduced_pt = bpt - Balls(center) val = sum( ZZ(j).factorial() * coeff(bpt) * polys[j](reduced_pt) for j, coeff in enumerate(post_transform)) return val
def __init__(self, point, dop=None, singular=None, **kwds): """ INPUT: - ``singular``: can be set to True to force this point to be considered a singular point, even if this cannot be checked (e.g. because we only have an enclosure) TESTS:: sage: from ore_algebra import * sage: from ore_algebra.analytic.path import Point sage: Dops, x, Dx = DifferentialOperators() sage: [Point(z, Dx) ....: for z in [1, 1/2, 1+I, QQbar(I), RIF(1/3), CIF(1/3), pi, ....: RDF(1), CDF(I), 0.5r, 0.5jr, 10r, QQbar(1), AA(1/3)]] [1, 1/2, I + 1, I, [0.333333333333333...], [0.333333333333333...], 3.141592653589794?, ~1.0000, ~1.0000*I, ~0.50000, ~0.50000*I, 10, 1, 1/3] sage: Point(sqrt(2), Dx).iv() [1.414...] sage: Point(RBF(0), (x-1)*x*Dx, singular=True).dist_to_sing() 1.000000000000000 """ SageObject.__init__(self) from sage.rings.complex_double import ComplexDoubleField_class from sage.rings.complex_field import ComplexField_class from sage.rings.complex_interval_field import ComplexIntervalField_class from sage.rings.real_double import RealDoubleField_class from sage.rings.real_mpfi import RealIntervalField_class from sage.rings.real_mpfr import RealField_class point = sage.structure.coerce.py_scalar_to_element(point) try: parent = point.parent() except AttributeError: raise TypeError("unexpected value for point: " + repr(point)) if isinstance(point, Point): self.value = point.value elif isinstance(parent, (RealBallField, ComplexBallField)): self.value = point elif isinstance(parent, number_field_base.NumberField): _, hom = good_number_field(point.parent()) self.value = hom(point) elif QQ.has_coerce_map_from(parent): self.value = QQ.coerce(point) elif QQbar.has_coerce_map_from(parent): alg = QQbar.coerce(point) NF, val, hom = alg.as_number_field_element() if NF is QQ: self.value = QQ.coerce(val) # parent may be ZZ else: embNF = number_field.NumberField(NF.polynomial(), NF.variable_name(), embedding=hom(NF.gen())) self.value = val.polynomial()(embNF.gen()) elif isinstance( parent, (RealField_class, RealDoubleField_class, RealIntervalField_class)): self.value = RealBallField(point.prec())(point) elif isinstance(parent, (ComplexField_class, ComplexDoubleField_class, ComplexIntervalField_class)): self.value = ComplexBallField(point.prec())(point) elif parent is sage.symbolic.ring.SR: try: return self.__init__(point.pyobject(), dop) except TypeError: pass try: return self.__init__(QQbar(point), dop) except (TypeError, ValueError, NotImplementedError): pass try: self.value = RLF(point) except (TypeError, ValueError): self.value = CLF(point) else: try: self.value = RLF.coerce(point) except TypeError: self.value = CLF.coerce(point) parent = self.value.parent() assert (isinstance( parent, (number_field_base.NumberField, RealBallField, ComplexBallField)) or parent is RLF or parent is CLF) if dop is None: # TBI if isinstance(point, Point): self.dop = point.dop else: self.dop = DifferentialOperator(dop.numerator()) self._force_singular = bool(singular) self.options = kwds