Regions in fundamental domains of period lattices

This module is used to represent sub-regions of a fundamental parallelogram of the period lattice of an elliptic curve, used in computing minimum height bounds.

In particular, these are the approximating sets S^{(v)} in section 3.2 of Thotsaphon Thongjunthug’s Ph.D. Thesis and paper [Tho2010].

AUTHORS:

  • Robert Bradshaw (2010): initial version

  • John Cremona (2014): added some docstrings and doctests

class sage.schemes.elliptic_curves.period_lattice_region.PeriodicRegion[source]

Bases: object

EXAMPLES:

sage: import numpy as np
sage: from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
sage: S = PeriodicRegion(CDF(2), CDF(2*I), np.zeros((4, 4)))
sage: S.plot()                                                              # needs sage.plot
Graphics object consisting of 1 graphics primitive
sage: data = np.zeros((4, 4))
sage: data[1,1] = True
sage: S = PeriodicRegion(CDF(2), CDF(2*I+1), data)
sage: S.plot()                                                              # needs sage.plot
Graphics object consisting of 5 graphics primitives
>>> from sage.all import *
>>> import numpy as np
>>> from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
>>> S = PeriodicRegion(CDF(Integer(2)), CDF(Integer(2)*I), np.zeros((Integer(4), Integer(4))))
>>> S.plot()                                                              # needs sage.plot
Graphics object consisting of 1 graphics primitive
>>> data = np.zeros((Integer(4), Integer(4)))
>>> data[Integer(1),Integer(1)] = True
>>> S = PeriodicRegion(CDF(Integer(2)), CDF(Integer(2)*I+Integer(1)), data)
>>> S.plot()                                                              # needs sage.plot
Graphics object consisting of 5 graphics primitives
border(raw=True)[source]

Return the boundary of this region as set of tile boundaries.

If raw is true, returns a list with respect to the internal bitmap, otherwise returns complex intervals covering the border.

EXAMPLES:

sage: import numpy as np
sage: from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
sage: data = np.zeros((4, 4))
sage: data[1, 1] = True
sage: PeriodicRegion(CDF(1), CDF(I), data).border()
[(1, 1, 0), (2, 1, 0), (1, 1, 1), (1, 2, 1)]
sage: PeriodicRegion(CDF(2), CDF(I-1/2), data).border()
[(1, 1, 0), (2, 1, 0), (1, 1, 1), (1, 2, 1)]

sage: PeriodicRegion(CDF(1), CDF(I), data).border(raw=False)
[0.25000000000000000? + 1.?*I,
 0.50000000000000000? + 1.?*I,
 1.? + 0.25000000000000000?*I,
 1.? + 0.50000000000000000?*I]
sage: PeriodicRegion(CDF(2), CDF(I-1/2), data).border(raw=False)
[0.3? + 1.?*I,
 0.8? + 1.?*I,
 1.? + 0.25000000000000000?*I,
 1.? + 0.50000000000000000?*I]

sage: data[1:3, 2] = True
sage: PeriodicRegion(CDF(1), CDF(I), data).border()
[(1, 1, 0), (2, 1, 0), (1, 1, 1), (1, 2, 0), (1, 3, 1), (3, 2, 0), (2, 2, 1), (2, 3, 1)]
>>> from sage.all import *
>>> import numpy as np
>>> from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
>>> data = np.zeros((Integer(4), Integer(4)))
>>> data[Integer(1), Integer(1)] = True
>>> PeriodicRegion(CDF(Integer(1)), CDF(I), data).border()
[(1, 1, 0), (2, 1, 0), (1, 1, 1), (1, 2, 1)]
>>> PeriodicRegion(CDF(Integer(2)), CDF(I-Integer(1)/Integer(2)), data).border()
[(1, 1, 0), (2, 1, 0), (1, 1, 1), (1, 2, 1)]

>>> PeriodicRegion(CDF(Integer(1)), CDF(I), data).border(raw=False)
[0.25000000000000000? + 1.?*I,
 0.50000000000000000? + 1.?*I,
 1.? + 0.25000000000000000?*I,
 1.? + 0.50000000000000000?*I]
>>> PeriodicRegion(CDF(Integer(2)), CDF(I-Integer(1)/Integer(2)), data).border(raw=False)
[0.3? + 1.?*I,
 0.8? + 1.?*I,
 1.? + 0.25000000000000000?*I,
 1.? + 0.50000000000000000?*I]

>>> data[Integer(1):Integer(3), Integer(2)] = True
>>> PeriodicRegion(CDF(Integer(1)), CDF(I), data).border()
[(1, 1, 0), (2, 1, 0), (1, 1, 1), (1, 2, 0), (1, 3, 1), (3, 2, 0), (2, 2, 1), (2, 3, 1)]
contract(corners=True)[source]

Opposite (but not inverse) of expand; removes neighbors of complement.

EXAMPLES:

sage: import numpy as np
sage: if int(np.version.short_version[0]) > 1:
....:     np.set_printoptions(legacy="1.25")
sage: from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
sage: data = np.zeros((10, 10))
sage: data[1:4,1:4] = True
sage: S = PeriodicRegion(CDF(1), CDF(I + 1/2), data)
sage: S.plot()                                                              # needs sage.plot
Graphics object consisting of 13 graphics primitives
sage: S.contract().plot()                                                   # needs sage.plot
Graphics object consisting of 5 graphics primitives
sage: S.contract().data.sum()
1
sage: S.contract().contract().is_empty()
True
>>> from sage.all import *
>>> import numpy as np
>>> if int(np.version.short_version[Integer(0)]) > Integer(1):
...     np.set_printoptions(legacy="1.25")
>>> from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
>>> data = np.zeros((Integer(10), Integer(10)))
>>> data[Integer(1):Integer(4),Integer(1):Integer(4)] = True
>>> S = PeriodicRegion(CDF(Integer(1)), CDF(I + Integer(1)/Integer(2)), data)
>>> S.plot()                                                              # needs sage.plot
Graphics object consisting of 13 graphics primitives
>>> S.contract().plot()                                                   # needs sage.plot
Graphics object consisting of 5 graphics primitives
>>> S.contract().data.sum()
1
>>> S.contract().contract().is_empty()
True
data[source]
ds()[source]

Return the sides of each parallelogram tile.

EXAMPLES:

sage: import numpy as np
sage: from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
sage: data = np.zeros((4, 4))
sage: S = PeriodicRegion(CDF(2), CDF(2*I), data, full=False)
sage: S.ds()
(0.5, 0.25*I)
sage: _ = S._ensure_full()
sage: S.ds()
(0.5, 0.25*I)

sage: data = np.zeros((8, 8))
sage: S = PeriodicRegion(CDF(1), CDF(I + 1/2), data)
sage: S.ds()
(0.125, 0.0625 + 0.125*I)
>>> from sage.all import *
>>> import numpy as np
>>> from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
>>> data = np.zeros((Integer(4), Integer(4)))
>>> S = PeriodicRegion(CDF(Integer(2)), CDF(Integer(2)*I), data, full=False)
>>> S.ds()
(0.5, 0.25*I)
>>> _ = S._ensure_full()
>>> S.ds()
(0.5, 0.25*I)

>>> data = np.zeros((Integer(8), Integer(8)))
>>> S = PeriodicRegion(CDF(Integer(1)), CDF(I + Integer(1)/Integer(2)), data)
>>> S.ds()
(0.125, 0.0625 + 0.125*I)
expand(corners=True)[source]

Return a region containing this region by adding all neighbors of internal tiles.

EXAMPLES:

sage: import numpy as np
sage: from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
sage: data = np.zeros((4, 4))
sage: data[1,1] = True
sage: S = PeriodicRegion(CDF(1), CDF(I + 1/2), data)
sage: S.plot()                                                              # needs sage.plot
Graphics object consisting of 5 graphics primitives
sage: S.expand().plot()                                                     # needs sage.plot
Graphics object consisting of 13 graphics primitives
sage: S.expand().data
array([[1, 1, 1, 0],
       [1, 1, 1, 0],
       [1, 1, 1, 0],
       [0, 0, 0, 0]], dtype=int8)
sage: S.expand(corners=False).plot()                                        # needs sage.plot
Graphics object consisting of 13 graphics primitives
sage: S.expand(corners=False).data
array([[0, 1, 0, 0],
       [1, 1, 1, 0],
       [0, 1, 0, 0],
       [0, 0, 0, 0]], dtype=int8)
>>> from sage.all import *
>>> import numpy as np
>>> from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
>>> data = np.zeros((Integer(4), Integer(4)))
>>> data[Integer(1),Integer(1)] = True
>>> S = PeriodicRegion(CDF(Integer(1)), CDF(I + Integer(1)/Integer(2)), data)
>>> S.plot()                                                              # needs sage.plot
Graphics object consisting of 5 graphics primitives
>>> S.expand().plot()                                                     # needs sage.plot
Graphics object consisting of 13 graphics primitives
>>> S.expand().data
array([[1, 1, 1, 0],
       [1, 1, 1, 0],
       [1, 1, 1, 0],
       [0, 0, 0, 0]], dtype=int8)
>>> S.expand(corners=False).plot()                                        # needs sage.plot
Graphics object consisting of 13 graphics primitives
>>> S.expand(corners=False).data
array([[0, 1, 0, 0],
       [1, 1, 1, 0],
       [0, 1, 0, 0],
       [0, 0, 0, 0]], dtype=int8)
full[source]
innermost_point()[source]

Return a point well inside the region, specifically the center of (one of) the last tile(s) to be removed on contraction.

EXAMPLES:

sage: import numpy as np
sage: from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
sage: data = np.zeros((10, 10))
sage: data[1:4, 1:4] = True
sage: data[1, 0:8] = True
sage: S = PeriodicRegion(CDF(1), CDF(I+1/2), data)
sage: S.innermost_point()
0.375 + 0.25*I
sage: S.plot() + point(S.innermost_point())                                 # needs sage.plot
Graphics object consisting of 24 graphics primitives
>>> from sage.all import *
>>> import numpy as np
>>> from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
>>> data = np.zeros((Integer(10), Integer(10)))
>>> data[Integer(1):Integer(4), Integer(1):Integer(4)] = True
>>> data[Integer(1), Integer(0):Integer(8)] = True
>>> S = PeriodicRegion(CDF(Integer(1)), CDF(I+Integer(1)/Integer(2)), data)
>>> S.innermost_point()
0.375 + 0.25*I
>>> S.plot() + point(S.innermost_point())                                 # needs sage.plot
Graphics object consisting of 24 graphics primitives
is_empty()[source]

Return whether this region is empty.

EXAMPLES:

sage: import numpy as np
sage: if int(np.version.short_version[0]) > 1:
....:     np.set_printoptions(legacy="1.25")
sage: from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
sage: data = np.zeros((4, 4))
sage: PeriodicRegion(CDF(2), CDF(2*I), data).is_empty()
True
sage: data[1,1] = True
sage: PeriodicRegion(CDF(2), CDF(2*I), data).is_empty()
False
>>> from sage.all import *
>>> import numpy as np
>>> if int(np.version.short_version[Integer(0)]) > Integer(1):
...     np.set_printoptions(legacy="1.25")
>>> from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
>>> data = np.zeros((Integer(4), Integer(4)))
>>> PeriodicRegion(CDF(Integer(2)), CDF(Integer(2)*I), data).is_empty()
True
>>> data[Integer(1),Integer(1)] = True
>>> PeriodicRegion(CDF(Integer(2)), CDF(Integer(2)*I), data).is_empty()
False
plot(**kwds)[source]

Plot this region in the fundamental lattice. If full is False, plots only the lower half. Note that the true nature of this region is periodic.

EXAMPLES:

sage: import numpy as np
sage: from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
sage: data = np.zeros((10, 10))
sage: data[2, 2:8] = True
sage: data[2:5, 2] = True
sage: data[3, 3] = True
sage: S = PeriodicRegion(CDF(1), CDF(I + 1/2), data)
sage: plot(S) + plot(S.expand(), rgbcolor=(1, 0, 1), thickness=2)           # needs sage.plot
Graphics object consisting of 46 graphics primitives
>>> from sage.all import *
>>> import numpy as np
>>> from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
>>> data = np.zeros((Integer(10), Integer(10)))
>>> data[Integer(2), Integer(2):Integer(8)] = True
>>> data[Integer(2):Integer(5), Integer(2)] = True
>>> data[Integer(3), Integer(3)] = True
>>> S = PeriodicRegion(CDF(Integer(1)), CDF(I + Integer(1)/Integer(2)), data)
>>> plot(S) + plot(S.expand(), rgbcolor=(Integer(1), Integer(0), Integer(1)), thickness=Integer(2))           # needs sage.plot
Graphics object consisting of 46 graphics primitives
refine(condition=None, times=1)[source]

Recursive function to refine the current tiling.

INPUT:

  • condition – function (default: None); if not None, only keep tiles in the refinement which satisfy the condition

  • times – integer (default: 1); the number of times to refine. Each refinement step halves the mesh size.

OUTPUT: the refined PeriodicRegion

EXAMPLES:

sage: import numpy as np
sage: from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
sage: data = np.zeros((4, 4))
sage: S = PeriodicRegion(CDF(2), CDF(2*I), data, full=False)
sage: S.ds()
(0.5, 0.25*I)
sage: S = S.refine()
sage: S.ds()
(0.25, 0.125*I)
sage: S = S.refine(2)
sage: S.ds()
(0.125, 0.0625*I)
>>> from sage.all import *
>>> import numpy as np
>>> from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
>>> data = np.zeros((Integer(4), Integer(4)))
>>> S = PeriodicRegion(CDF(Integer(2)), CDF(Integer(2)*I), data, full=False)
>>> S.ds()
(0.5, 0.25*I)
>>> S = S.refine()
>>> S.ds()
(0.25, 0.125*I)
>>> S = S.refine(Integer(2))
>>> S.ds()
(0.125, 0.0625*I)
verify(condition)[source]

Given a condition that should hold for every line segment on the boundary, verify that it actually does so.

INPUT:

  • condition – boolean-valued function on \(\CC\)

OUTPUT:

boolean according to whether the condition holds for all lines on the boundary.

EXAMPLES:

sage: import numpy as np
sage: from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
sage: data = np.zeros((4, 4))
sage: data[1, 1] = True
sage: S = PeriodicRegion(CDF(1), CDF(I), data)
sage: S.border()
[(1, 1, 0), (2, 1, 0), (1, 1, 1), (1, 2, 1)]
sage: condition = lambda z: z.real().abs()<1/2
sage: S.verify(condition)
False
sage: condition = lambda z: z.real().abs()<1
sage: S.verify(condition)
True
>>> from sage.all import *
>>> import numpy as np
>>> from sage.schemes.elliptic_curves.period_lattice_region import PeriodicRegion
>>> data = np.zeros((Integer(4), Integer(4)))
>>> data[Integer(1), Integer(1)] = True
>>> S = PeriodicRegion(CDF(Integer(1)), CDF(I), data)
>>> S.border()
[(1, 1, 0), (2, 1, 0), (1, 1, 1), (1, 2, 1)]
>>> condition = lambda z: z.real().abs()<Integer(1)/Integer(2)
>>> S.verify(condition)
False
>>> condition = lambda z: z.real().abs()<Integer(1)
>>> S.verify(condition)
True
w1[source]
w2[source]