Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagesmc
Path: blob/master/src/sage/schemes/toric/variety.py
8820 views
1
# -*- coding: utf-8 -*-
2
r"""
3
Toric varieties
4
5
This module provides support for (normal) toric varieties, corresponding to
6
:class:`rational polyhedral fans <sage.geometry.fan.RationalPolyhedralFan>`.
7
See also :mod:`~sage.schemes.toric.fano_variety` for a more
8
restrictive class of (weak) Fano toric varieties.
9
10
An **excellent reference on toric varieties** is the book "Toric
11
Varieties" by David A. Cox, John B. Little, and Hal Schenck
12
[CLS]_.
13
14
The interface to this module is provided through functions
15
:func:`AffineToricVariety` and :func:`ToricVariety`, although you may
16
also be interested in :func:`normalize_names`.
17
18
.. NOTE::
19
20
We do NOT build "general toric varieties" from affine toric varieties.
21
Instead, we are using the quotient representation of toric varieties with
22
the homogeneous coordinate ring (a.k.a. Cox's ring or the total coordinate
23
ring). This description works best for simplicial fans of the full
24
dimension.
25
26
REFERENCES:
27
28
.. [CLS]
29
David A. Cox, John B. Little, Hal Schenck,
30
"Toric Varieties", Graduate Studies in Mathematics,
31
Amer. Math. Soc., Providence, RI, 2011
32
33
AUTHORS:
34
35
- Andrey Novoseltsev (2010-05-17): initial version.
36
37
- Volker Braun (2010-07-24): Cohomology and characteristic classes added.
38
39
EXAMPLES:
40
41
We start with constructing the affine plane as an affine toric variety. First,
42
we need to have a corresponding cone::
43
44
sage: quadrant = Cone([(1,0), (0,1)])
45
46
If you don't care about variable names and the base field, that's all we need
47
for now::
48
49
sage: A2 = AffineToricVariety(quadrant)
50
sage: A2
51
2-d affine toric variety
52
sage: origin = A2(0,0)
53
sage: origin
54
[0 : 0]
55
56
Only affine toric varieties have points whose (homogeneous) coordinates
57
are all zero. ::
58
59
sage: parent(origin)
60
Set of rational points of 2-d affine toric variety
61
62
As you can see, by default toric varieties live over the field of rational
63
numbers::
64
65
sage: A2.base_ring()
66
Rational Field
67
68
While usually toric varieties are considered over the field of complex
69
numbers, for computational purposes it is more convenient to work with fields
70
that have exact representation on computers. You can also always do ::
71
72
sage: C2 = AffineToricVariety(quadrant, base_field=CC)
73
sage: C2.base_ring()
74
Complex Field with 53 bits of precision
75
sage: C2(1,2+i)
76
[1.00000000000000 : 2.00000000000000 + 1.00000000000000*I]
77
78
or even ::
79
80
sage: F = CC["a, b"].fraction_field()
81
sage: F.inject_variables()
82
Defining a, b
83
sage: A2 = AffineToricVariety(quadrant, base_field=F)
84
sage: A2(a,b)
85
[a : b]
86
87
OK, if you need to work only with affine spaces,
88
:func:`~sage.schemes.affine.affine_space.AffineSpace` may be a better way to
89
construct them. Our next example is the product of two projective lines
90
realized as the toric variety associated to the
91
:func:`face fan <sage.geometry.fan.FaceFan>` of the "diamond"::
92
93
sage: diamond = lattice_polytope.octahedron(2)
94
sage: diamond.vertices()
95
[ 1 0 -1 0]
96
[ 0 1 0 -1]
97
sage: fan = FaceFan(diamond)
98
sage: P1xP1 = ToricVariety(fan)
99
sage: P1xP1
100
2-d toric variety covered by 4 affine patches
101
sage: P1xP1.fan().rays()
102
N( 1, 0),
103
N( 0, 1),
104
N(-1, 0),
105
N( 0, -1)
106
in 2-d lattice N
107
sage: P1xP1.gens()
108
(z0, z1, z2, z3)
109
110
We got four coordinates - two for each of the projective lines, but their
111
names are perhaps not very well chosen. Let's make `(x,y)` to be coordinates
112
on the first line and `(s,t)` on the second one::
113
114
sage: P1xP1 = ToricVariety(fan, coordinate_names="x s y t")
115
sage: P1xP1.gens()
116
(x, s, y, t)
117
118
Now, if we want to define subschemes of this variety, the defining polynomials
119
must be homogeneous in each of these pairs::
120
121
sage: P1xP1.inject_variables()
122
Defining x, s, y, t
123
sage: P1xP1.subscheme(x)
124
Closed subscheme of 2-d toric variety
125
covered by 4 affine patches defined by:
126
x
127
sage: P1xP1.subscheme(x^2 + y^2)
128
Closed subscheme of 2-d toric variety
129
covered by 4 affine patches defined by:
130
x^2 + y^2
131
sage: P1xP1.subscheme(x^2 + s^2)
132
Traceback (most recent call last):
133
...
134
ValueError: x^2 + s^2 is not homogeneous
135
on 2-d toric variety covered by 4 affine patches!
136
sage: P1xP1.subscheme([x^2*s^2 + x*y*t^2 +y^2*t^2, s^3 + t^3])
137
Closed subscheme of 2-d toric variety
138
covered by 4 affine patches defined by:
139
x^2*s^2 + x*y*t^2 + y^2*t^2,
140
s^3 + t^3
141
142
While we don't build toric varieties from affine toric varieties, we still can
143
access the "building pieces"::
144
145
sage: patch = P1xP1.affine_patch(2)
146
sage: patch
147
2-d affine toric variety
148
sage: patch.fan().rays()
149
N(1, 0),
150
N(0, 1)
151
in 2-d lattice N
152
sage: patch.embedding_morphism()
153
Scheme morphism:
154
From: 2-d affine toric variety
155
To: 2-d toric variety covered by 4 affine patches
156
Defn: Defined on coordinates by sending [x : s] to
157
[x : s : 1 : 1]
158
159
The patch above was specifically chosen to coincide with our representation of
160
the affine plane before, but you can get the other three patches as well.
161
(While any cone of a fan will correspond to an affine toric variety, the main
162
interest is usually in the generating fans as "the biggest" affine
163
subvarieties, and these are precisely the patches that you can get from
164
:meth:`~ToricVariety_field.affine_patch`.)
165
166
All two-dimensional toric varieties are "quite nice" because any
167
two-dimensional cone is generated by exactly two rays. From the point of view
168
of the corresponding toric varieties, this means that they have at worst
169
quotient singularities::
170
171
sage: P1xP1.is_orbifold()
172
True
173
sage: P1xP1.is_smooth()
174
True
175
sage: TV = ToricVariety(NormalFan(diamond))
176
sage: TV.fan().rays()
177
N(-1, 1),
178
N( 1, 1),
179
N(-1, -1),
180
N( 1, -1)
181
in 2-d lattice N
182
sage: TV.is_orbifold()
183
True
184
sage: TV.is_smooth()
185
False
186
187
In higher dimensions worse things can happen::
188
189
sage: TV3 = ToricVariety(NormalFan(lattice_polytope.octahedron(3)))
190
sage: TV3.fan().rays()
191
N(-1, -1, 1),
192
N( 1, -1, 1),
193
N(-1, 1, 1),
194
N( 1, 1, 1),
195
N(-1, -1, -1),
196
N( 1, -1, -1),
197
N(-1, 1, -1),
198
N( 1, 1, -1)
199
in 3-d lattice N
200
sage: TV3.is_orbifold()
201
False
202
203
Fortunately, we can perform a (partial) resolution::
204
205
sage: TV3_res = TV3.resolve_to_orbifold()
206
sage: TV3_res.is_orbifold()
207
True
208
sage: TV3_res.fan().ngenerating_cones()
209
12
210
sage: TV3.fan().ngenerating_cones()
211
6
212
213
In this example we had to double the number of affine patches. The result is
214
still singular::
215
216
sage: TV3_res.is_smooth()
217
False
218
219
You can resolve it further using :meth:`~ToricVariety_field.resolve` method,
220
but (at least for now) you will have to specify which rays should be inserted
221
into the fan. See also
222
:func:`~sage.schemes.toric.fano_variety.CPRFanoToricVariety`,
223
which can construct some other "nice partial resolutions."
224
225
The intersection theory on toric varieties is very well understood,
226
and there are explicit algorithms to compute many quantities of
227
interest. The most important tools are the :class:`cohomology ring
228
<CohomologyRing>` and the :mod:`Chow group
229
<sage.schemes.toric.chow_group>`. For `d`-dimensional compact
230
toric varieties with at most orbifold singularities, the rational
231
cohomology ring `H^*(X,\QQ)` and the rational Chow ring `A^*(X,\QQ) =
232
A_{d-*}(X)\otimes \QQ` are isomorphic except for a doubling in
233
degree. More precisely, the Chow group has the same rank
234
235
.. math::
236
237
A_{d-k}(X) \otimes \QQ \simeq H^{2k}(X,\QQ)
238
239
and the intersection in of Chow cycles matches the cup product in
240
cohomology.
241
242
In this case, you should work with the cohomology ring description
243
because it is much faster. For example, here is a weighted projective
244
space with a curve of `\ZZ_3`-orbifold singularities::
245
246
sage: P4_11133 = toric_varieties.P4_11133()
247
sage: P4_11133.is_smooth(), P4_11133.is_orbifold()
248
(False, True)
249
sage: cone = P4_11133.fan(3)[8]
250
sage: cone.is_smooth(), cone.is_simplicial()
251
(False, True)
252
sage: HH = P4_11133.cohomology_ring(); HH
253
Rational cohomology ring of a 4-d CPR-Fano toric variety covered by 5 affine patches
254
sage: P4_11133.cohomology_basis()
255
(([1],), ([z4],), ([z4^2],), ([z4^3],), ([z4^4],))
256
257
Every cone defines a torus orbit closure, and hence a (co)homology class::
258
259
sage: HH.gens()
260
([3*z4], [3*z4], [z4], [z4], [z4])
261
sage: map(HH, P4_11133.fan(1))
262
[[3*z4], [3*z4], [z4], [z4], [z4]]
263
sage: map(HH, P4_11133.fan(4) )
264
[[9*z4^4], [9*z4^4], [9*z4^4], [9*z4^4], [9*z4^4]]
265
sage: HH(cone)
266
[3*z4^3]
267
268
We can compute intersection numbers by integrating top-dimensional
269
cohomology classes::
270
271
sage: D = P4_11133.divisor(0)
272
sage: HH(D)
273
[3*z4]
274
sage: P4_11133.integrate( HH(D)^4 )
275
9
276
sage: P4_11133.integrate( HH(D) * HH(cone) )
277
1
278
279
Although computationally less efficient, we can do the same
280
computations with the rational Chow group::
281
282
sage: AA = P4_11133.Chow_group(QQ)
283
sage: map(AA, P4_11133.fan(1)) # long time (5s on sage.math, 2012)
284
[( 0 | 0 | 0 | 3 | 0 ), ( 0 | 0 | 0 | 3 | 0 ), ( 0 | 0 | 0 | 1 | 0 ), ( 0 | 0 | 0 | 1 | 0 ), ( 0 | 0 | 0 | 1 | 0 )]
285
sage: map(AA, P4_11133.fan(4)) # long time (5s on sage.math, 2012)
286
[( 1 | 0 | 0 | 0 | 0 ), ( 1 | 0 | 0 | 0 | 0 ), ( 1 | 0 | 0 | 0 | 0 ), ( 1 | 0 | 0 | 0 | 0 ), ( 1 | 0 | 0 | 0 | 0 )]
287
sage: AA(cone).intersection_with_divisor(D) # long time (4s on sage.math, 2013)
288
( 1 | 0 | 0 | 0 | 0 )
289
sage: AA(cone).intersection_with_divisor(D).count_points() # long time
290
1
291
292
The real advantage of the Chow group is that
293
294
* it works just as well over `\ZZ`, so torsion information is also
295
easily available, and
296
297
* its combinatorial description also works over worse-than-orbifold
298
singularities. By contrast, the cohomology groups can become very
299
complicated to compute in this case, and one usually only has a
300
spectral sequence but no toric algorithm.
301
302
Below you will find detailed descriptions of available functions. If you are
303
familiar with toric geometry, you will likely see that many important objects
304
and operations are unavailable. However, this module is under active
305
development and hopefully will improve in future releases of Sage. If there
306
are some particular features that you would like to see implemented ASAP,
307
please consider reporting them to the Sage Development Team or even
308
implementing them on your own as a patch for inclusion!
309
"""
310
311
312
#*****************************************************************************
313
# Copyright (C) 2010 Volker Braun <[email protected]>
314
# Copyright (C) 2010 Andrey Novoseltsev <[email protected]>
315
# Copyright (C) 2010 William Stein <[email protected]>
316
#
317
# Distributed under the terms of the GNU General Public License (GPL)
318
# as published by the Free Software Foundation; either version 2 of
319
# the License, or (at your option) any later version.
320
# http://www.gnu.org/licenses/
321
#*****************************************************************************
322
323
import sys
324
325
from sage.functions.all import factorial
326
from sage.geometry.cone import Cone, is_Cone
327
from sage.geometry.fan import Fan
328
from sage.matrix.all import matrix
329
from sage.misc.all import latex, prod, uniq, cached_method
330
from sage.structure.unique_representation import UniqueRepresentation
331
from sage.modules.free_module_element import vector
332
from sage.rings.all import Infinity, PolynomialRing, ZZ, QQ
333
from sage.rings.quotient_ring_element import QuotientRingElement
334
from sage.rings.quotient_ring import QuotientRing_generic
335
from sage.schemes.affine.affine_space import AffineSpace
336
from sage.schemes.generic.ambient_space import AmbientSpace
337
from sage.schemes.toric.homset import SchemeHomset_points_toric_field
338
from sage.categories.fields import Fields
339
from sage.misc.cachefunc import ClearCacheOnPickle
340
_Fields = Fields()
341
342
343
# Default prefix for indexed coordinates
344
DEFAULT_PREFIX = "z"
345
346
347
def is_ToricVariety(x):
348
r"""
349
Check if ``x`` is a toric variety.
350
351
INPUT:
352
353
- ``x`` -- anything.
354
355
OUTPUT:
356
357
- ``True`` if ``x`` is a :class:`toric variety <ToricVariety_field>` and
358
``False`` otherwise.
359
360
.. NOTE::
361
362
While projective spaces are toric varieties mathematically, they are
363
not toric varieties in Sage due to efficiency considerations, so this
364
function will return ``False``.
365
366
EXAMPLES::
367
368
sage: from sage.schemes.toric.variety import is_ToricVariety
369
sage: is_ToricVariety(1)
370
False
371
sage: fan = FaceFan(lattice_polytope.octahedron(2))
372
sage: P = ToricVariety(fan)
373
sage: P
374
2-d toric variety covered by 4 affine patches
375
sage: is_ToricVariety(P)
376
True
377
sage: is_ToricVariety(ProjectiveSpace(2))
378
False
379
"""
380
return isinstance(x, ToricVariety_field)
381
382
383
def ToricVariety(fan,
384
coordinate_names=None,
385
names=None,
386
coordinate_indices=None,
387
base_ring=QQ, base_field=None):
388
r"""
389
Construct a toric variety.
390
391
INPUT:
392
393
- ``fan`` -- :class:`rational polyhedral fan
394
<sage.geometry.fan.RationalPolyhedralFan>`;
395
396
- ``coordinate_names`` -- names of variables for the coordinate ring, see
397
:func:`normalize_names` for acceptable formats. If not given, indexed
398
variable names will be created automatically;
399
400
- ``names`` -- an alias of ``coordinate_names`` for internal
401
use. You may specify either ``names`` or ``coordinate_names``,
402
but not both;
403
404
- ``coordinate_indices`` -- list of integers, indices for indexed
405
variables. If not given, the index of each variable will coincide with
406
the index of the corresponding ray of the fan;
407
408
- ``base_ring`` -- base ring of the toric variety (default:
409
`\QQ`). Must be a field.
410
411
- ``base_field`` -- alias for ``base_ring``. Takes precedence if
412
both are specified.
413
414
OUTPUT:
415
416
- :class:`toric variety <ToricVariety_field>`.
417
418
EXAMPLES:
419
420
We will create the product of two projective lines::
421
422
sage: fan = FaceFan(lattice_polytope.octahedron(2))
423
sage: fan.rays()
424
N( 1, 0),
425
N( 0, 1),
426
N(-1, 0),
427
N( 0, -1)
428
in 2-d lattice N
429
sage: P1xP1 = ToricVariety(fan)
430
sage: P1xP1.gens()
431
(z0, z1, z2, z3)
432
433
Let's create some points::
434
435
sage: P1xP1(1,1,1,1)
436
[1 : 1 : 1 : 1]
437
sage: P1xP1(0,1,1,1)
438
[0 : 1 : 1 : 1]
439
sage: P1xP1(0,1,0,1)
440
Traceback (most recent call last):
441
...
442
TypeError: coordinates (0, 1, 0, 1)
443
are in the exceptional set!
444
445
We cannot set to zero both coordinates of the same projective line!
446
447
Let's change the names of the variables. We have to re-create our toric
448
variety::
449
450
sage: P1xP1 = ToricVariety(fan, "x s y t")
451
sage: P1xP1.gens()
452
(x, s, y, t)
453
454
Now `(x, y)` correspond to one line and `(s, t)` to the other one. ::
455
456
sage: P1xP1.inject_variables()
457
Defining x, s, y, t
458
sage: P1xP1.subscheme(x*s-y*t)
459
Closed subscheme of 2-d toric variety
460
covered by 4 affine patches defined by:
461
x*s - y*t
462
463
Here is a shorthand for defining the toric variety and homogeneous
464
coordinates in one go::
465
466
sage: P1xP1.<a,b,c,d> = ToricVariety(fan)
467
sage: (a^2+b^2) * (c+d)
468
a^2*c + b^2*c + a^2*d + b^2*d
469
"""
470
if base_field is not None:
471
base_ring = base_field
472
if names is not None:
473
if coordinate_names is not None:
474
raise ValueError('You must not specify both coordinate_names and names!')
475
coordinate_names = names
476
if base_ring not in _Fields:
477
raise TypeError("need a field to construct a toric variety!\n Got %s"
478
% base_ring)
479
return ToricVariety_field(fan, coordinate_names, coordinate_indices,
480
base_ring)
481
482
483
def AffineToricVariety(cone, *args, **kwds):
484
r"""
485
Construct an affine toric variety.
486
487
INPUT:
488
489
- ``cone`` -- :class:`strictly convex rational polyhedral cone
490
<sage.geometry.cone.ConvexRationalPolyhedralCone>`.
491
492
This cone will be used to construct a :class:`rational polyhedral fan
493
<sage.geometry.fan.RationalPolyhedralFan>`, which will be passed to
494
:func:`ToricVariety` with the rest of positional and keyword arguments.
495
496
OUTPUT:
497
498
- :class:`toric variety <ToricVariety_field>`.
499
500
.. NOTE::
501
502
The generating rays of the fan of this variety are guaranteed to be
503
listed in the same order as the rays of the original cone.
504
505
EXAMPLES:
506
507
We will create the affine plane as an affine toric variety::
508
509
sage: quadrant = Cone([(1,0), (0,1)])
510
sage: A2 = AffineToricVariety(quadrant)
511
sage: origin = A2(0,0)
512
sage: origin
513
[0 : 0]
514
sage: parent(origin)
515
Set of rational points of 2-d affine toric variety
516
517
Only affine toric varieties have points whose (homogeneous) coordinates
518
are all zero.
519
"""
520
if not cone.is_strictly_convex():
521
raise ValueError("affine toric varieties are defined for strictly "
522
"convex cones only!")
523
# We make sure that Fan constructor does not meddle with the order of
524
# rays, this is very important for affine patches construction
525
fan = Fan([tuple(range(cone.nrays()))], cone.rays(),
526
check=False, normalize=False)
527
return ToricVariety(fan, *args, **kwds)
528
529
530
class ToricVariety_field(ClearCacheOnPickle, AmbientSpace):
531
r"""
532
Construct a toric variety associated to a rational polyhedral fan.
533
534
.. WARNING::
535
536
This class does not perform any checks of correctness of input. Use
537
:func:`ToricVariety` and :func:`AffineToricVariety` to construct toric
538
varieties.
539
540
INPUT:
541
542
- ``fan`` -- :class:`rational polyhedral fan
543
<sage.geometry.fan.RationalPolyhedralFan>`;
544
545
- ``coordinate_names`` -- names of variables, see :func:`normalize_names`
546
for acceptable formats. If ``None``, indexed variable names will be
547
created automatically;
548
549
- ``coordinate_indices`` -- list of integers, indices for indexed
550
variables. If ``None``, the index of each variable will coincide with
551
the index of the corresponding ray of the fan;
552
553
- ``base_field`` -- base field of the toric variety.
554
555
OUTPUT:
556
557
- :class:`toric variety <ToricVariety_field>`.
558
559
TESTS::
560
561
sage: fan = FaceFan(lattice_polytope.octahedron(2))
562
sage: P1xP1 = ToricVariety(fan)
563
"""
564
565
def __init__(self, fan, coordinate_names, coordinate_indices, base_field):
566
r"""
567
See :class:`ToricVariety_field` for documentation.
568
569
TESTS::
570
571
sage: fan = FaceFan(lattice_polytope.octahedron(2))
572
sage: P1xP1 = ToricVariety(fan)
573
"""
574
self._fan = fan
575
super(ToricVariety_field, self).__init__(fan.lattice_dim(),
576
base_field)
577
self._torus_factor_dim = fan.lattice_dim() - fan.dim()
578
coordinate_names = normalize_names(coordinate_names,
579
fan.nrays() + self._torus_factor_dim, DEFAULT_PREFIX,
580
coordinate_indices, return_prefix=True)
581
# Save the prefix for use in resolutions
582
self._coordinate_prefix = coordinate_names.pop()
583
self._assign_names(names=coordinate_names, normalize=False)
584
585
def __cmp__(self, right):
586
r"""
587
Compare ``self`` and ``right``.
588
589
INPUT:
590
591
- ``right`` -- anything.
592
593
OUTPUT:
594
595
- 0 if ``right`` is of the same type as ``self``, their fans are the
596
same, names of variables are the same and stored in the same order,
597
and base fields are the same. 1 or -1 otherwise.
598
599
TESTS::
600
601
sage: fan = FaceFan(lattice_polytope.octahedron(2))
602
sage: P1xP1 = ToricVariety(fan)
603
sage: P1xP1a = ToricVariety(fan, "x s y t")
604
sage: P1xP1b = ToricVariety(fan)
605
sage: cmp(P1xP1, P1xP1a)
606
1
607
sage: cmp(P1xP1a, P1xP1)
608
-1
609
sage: cmp(P1xP1, P1xP1b)
610
0
611
sage: P1xP1 is P1xP1b
612
False
613
sage: cmp(P1xP1, 1) * cmp(1, P1xP1)
614
-1
615
"""
616
c = cmp(type(self), type(right))
617
if c:
618
return c
619
return cmp([self.fan(),
620
self.variable_names(),
621
self.base_ring()],
622
[right.fan(),
623
right.variable_names(),
624
right.base_ring()])
625
626
def _an_element_(self):
627
r"""
628
Construct an element of ``self``.
629
630
This function is needed (in particular) for the test framework.
631
632
OUTPUT:
633
634
- a point of ``self`` with coordinates [1 : 2: ... : n].
635
636
TESTS::
637
638
sage: fan = FaceFan(lattice_polytope.octahedron(2))
639
sage: P1xP1 = ToricVariety(fan)
640
sage: P1xP1._an_element_()
641
[1 : 2 : 3 : 4]
642
"""
643
return self(range(1, self.ngens() + 1))
644
645
def _check_satisfies_equations(self, coordinates):
646
r"""
647
Check if ``coordinates`` define a valid point of ``self``.
648
649
INPUT:
650
651
- ``coordinates`` -- list of elements of the base field of ``self``.
652
653
OUTPUT:
654
655
- ``True`` if ``coordinates`` do define a valid point of ``self``,
656
otherwise a ``TypeError`` or ``ValueError`` exception is raised.
657
658
TESTS::
659
660
sage: fan = FaceFan(lattice_polytope.octahedron(2))
661
sage: P1xP1 = ToricVariety(fan)
662
sage: P1xP1._check_satisfies_equations([1,1,1,1])
663
True
664
sage: P1xP1._check_satisfies_equations([0,0,1,1])
665
True
666
sage: P1xP1._check_satisfies_equations([0,1,0,1])
667
Traceback (most recent call last):
668
...
669
TypeError: coordinates (0, 1, 0, 1)
670
are in the exceptional set!
671
sage: P1xP1._check_satisfies_equations([1,1,1])
672
Traceback (most recent call last):
673
...
674
TypeError: coordinates (1, 1, 1) must have 4 components!
675
sage: P1xP1._check_satisfies_equations(1)
676
Traceback (most recent call last):
677
...
678
TypeError: 1 can not be used as coordinates!
679
Use a list or a tuple.
680
sage: P1xP1._check_satisfies_equations([1,1,1,fan])
681
Traceback (most recent call last):
682
...
683
TypeError: coordinate Rational polyhedral fan
684
in 2-d lattice N is not an element of Rational Field!
685
"""
686
try:
687
coordinates = tuple(coordinates)
688
except TypeError:
689
raise TypeError("%s can not be used as coordinates! "
690
"Use a list or a tuple." % coordinates)
691
n = self.ngens()
692
if len(coordinates) != n:
693
raise TypeError("coordinates %s must have %d components!"
694
% (coordinates, n))
695
base_field = self.base_ring()
696
for coordinate in coordinates:
697
if coordinate not in base_field:
698
raise TypeError("coordinate %s is not an element of %s!"
699
% (coordinate, base_field))
700
zero_positions = set(position
701
for position, coordinate in enumerate(coordinates)
702
if coordinate == 0)
703
if not zero_positions:
704
return True
705
for i in range(n - self._torus_factor_dim, n):
706
if i in zero_positions:
707
raise ValueError("coordinates on the torus factor cannot be "
708
"zero! Got %s" % str(coordinates))
709
if len(zero_positions) == 1:
710
return True
711
fan = self.fan()
712
possible_charts = set(fan._ray_to_cones(zero_positions.pop()))
713
for i in zero_positions:
714
possible_charts.intersection_update(fan._ray_to_cones(i))
715
if possible_charts:
716
return True # All zeros are inside one generating cone
717
raise TypeError("coordinates %s are in the exceptional set!"
718
% str(coordinates)) # Need str, coordinates is a tuple
719
720
def _point_homset(self, *args, **kwds):
721
r"""
722
Construct a Hom-set for ``self``.
723
724
INPUT:
725
726
- same as for
727
:class:`~sage.schemes.generic.homset.SchemeHomset_points_toric_field`.
728
729
OUPUT:
730
731
-
732
:class:`~sage.schemes.generic.homset.SchemeHomset_points_toric_field`.
733
734
TESTS::
735
736
sage: fan = FaceFan(lattice_polytope.octahedron(2))
737
sage: P1xP1 = ToricVariety(fan)
738
sage: P1xP1._point_homset(Spec(QQ), P1xP1)
739
Set of rational points of 2-d toric variety
740
covered by 4 affine patches
741
"""
742
return SchemeHomset_points_toric_field(*args, **kwds)
743
744
def _latex_(self):
745
r"""
746
Return a LaTeX representation of ``self``.
747
748
OUTPUT:
749
750
- string.
751
752
TESTS::
753
754
sage: fan = FaceFan(lattice_polytope.octahedron(2))
755
sage: P1xP1 = ToricVariety(fan)
756
sage: P1xP1._latex_()
757
'\\mathbb{X}_{\\Sigma^{2}}'
758
"""
759
return r"\mathbb{X}_{%s}" % latex(self.fan())
760
761
def _latex_generic_point(self, coordinates=None):
762
"""
763
Return a LaTeX representation of a point of ``self``.
764
765
INPUT:
766
767
- ``coordinates`` -- list of coordinates of a point of ``self``.
768
If not given, names of coordinates of ``self`` will be used.
769
770
OUTPUT:
771
772
- string.
773
774
EXAMPLES::
775
776
sage: fan = FaceFan(lattice_polytope.octahedron(2))
777
sage: P1xP1 = ToricVariety(fan)
778
sage: P1xP1._latex_generic_point()
779
'\\left[z_{0} : z_{1} : z_{2} : z_{3}\\right]'
780
sage: P1xP1._latex_generic_point([1,2,3,4])
781
'\\left[1 : 2 : 3 : 4\\right]'
782
"""
783
if coordinates is None:
784
coordinates = self.gens()
785
return r"\left[%s\right]" % (" : ".join(str(latex(coord))
786
for coord in coordinates))
787
788
def _point(self, *args, **kwds):
789
r"""
790
Construct a point of ``self``.
791
792
INPUT:
793
794
- same as for
795
:class:`~sage.schemes.generic.morphism.SchemeMorphism_point_toric_field`.
796
797
OUPUT:
798
799
:class:`~sage.schemes.generic.morphism.SchemeMorphism_point_toric_field`.
800
801
TESTS::
802
803
sage: fan = FaceFan(lattice_polytope.octahedron(2))
804
sage: P1xP1 = ToricVariety(fan)
805
sage: P1xP1._point(P1xP1, [1,2,3,4])
806
[1 : 2 : 3 : 4]
807
"""
808
from sage.schemes.toric.morphism import SchemeMorphism_point_toric_field
809
return SchemeMorphism_point_toric_field(*args, **kwds)
810
811
def _homset(self, *args, **kwds):
812
r"""
813
Return the homset between two toric varieties.
814
815
INPUT:
816
817
Same as :class:`sage.schemes.generic.homset.SchemeHomset_generic`.
818
819
OUTPUT:
820
821
A :class:`sage.schemes.toric.homset.SchemeHomset_toric_variety`.
822
823
EXAMPLES::
824
825
sage: P1xP1 = toric_varieties.P1xP1()
826
sage: P1 = toric_varieties.P1()
827
sage: hom_set = P1xP1.Hom(P1); hom_set
828
Set of morphisms
829
From: 2-d CPR-Fano toric variety covered by 4 affine patches
830
To: 1-d CPR-Fano toric variety covered by 2 affine patches
831
sage: type(hom_set)
832
<class 'sage.schemes.toric.homset.SchemeHomset_toric_variety_with_category'>
833
834
This is also the Hom-set for algebraic subschemes of toric varieties::
835
836
sage: P1xP1.inject_variables()
837
Defining s, t, x, y
838
sage: P1 = P1xP1.subscheme(s-t)
839
sage: hom_set = P1xP1.Hom(P1)
840
sage: hom_set([s,s,x,y])
841
Scheme morphism:
842
From: 2-d CPR-Fano toric variety covered by 4 affine patches
843
To: Closed subscheme of 2-d CPR-Fano toric variety covered by 4 affine patches defined by:
844
s - t
845
Defn: Defined on coordinates by sending [s : t : x : y] to
846
[s : s : x : y]
847
848
sage: hom_set = P1.Hom(P1)
849
sage: hom_set([s,s,x,y])
850
Scheme endomorphism of Closed subscheme of 2-d CPR-Fano toric
851
variety covered by 4 affine patches defined by:
852
s - t
853
Defn: Defined on coordinates by sending [s : t : x : y] to
854
[t : t : x : y]
855
"""
856
from sage.schemes.toric.homset import SchemeHomset_toric_variety
857
return SchemeHomset_toric_variety(*args, **kwds)
858
859
def _repr_(self):
860
r"""
861
Return a string representation of ``self``.
862
863
OUTPUT:
864
865
- string.
866
867
TESTS::
868
869
sage: fan = FaceFan(lattice_polytope.octahedron(2))
870
sage: P1xP1 = ToricVariety(fan)
871
sage: P1xP1._repr_()
872
'2-d toric variety covered by 4 affine patches'
873
"""
874
result = "%d-d" % self.dimension_relative()
875
if self.fan().ngenerating_cones() == 1:
876
result += " affine toric variety"
877
else:
878
result += (" toric variety covered by %d affine patches"
879
% self.fan().ngenerating_cones())
880
return result
881
882
def _repr_generic_point(self, coordinates=None):
883
r"""
884
Return a string representation of a point of ``self``.
885
886
INPUT:
887
888
- ``coordinates`` -- list of coordinates of a point of ``self``.
889
If not given, names of coordinates of ``self`` will be used.
890
891
OUTPUT:
892
893
- string.
894
895
EXAMPLES::
896
897
sage: fan = FaceFan(lattice_polytope.octahedron(2))
898
sage: P1xP1 = ToricVariety(fan)
899
sage: P1xP1._repr_generic_point()
900
'[z0 : z1 : z2 : z3]'
901
sage: P1xP1._repr_generic_point([1,2,3,4])
902
'[1 : 2 : 3 : 4]'
903
"""
904
if coordinates is None:
905
coordinates = self.gens()
906
return "[%s]" % (" : ".join(str(coord) for coord in coordinates))
907
908
def _validate(self, polynomials):
909
"""
910
Check if ``polynomials`` define valid functions on ``self``.
911
912
Since this is a toric variety, polynomials must be homogeneous in the
913
total coordinate ring of ``self``.
914
915
INPUT:
916
917
- ``polynomials`` -- list of polynomials in the coordinate ring of
918
``self`` (this function does not perform any conversions).
919
920
OUTPUT:
921
922
- ``polynomials`` (the input parameter without any modifications) if
923
``polynomials`` do define valid polynomial functions on ``self``,
924
otherwise a ``ValueError`` exception is raised.
925
926
TESTS:
927
928
We will use the product of two projective lines with coordinates
929
`(x, y)` for one and `(s, t)` for the other::
930
931
sage: fan = FaceFan(lattice_polytope.octahedron(2))
932
sage: fan.rays()
933
N( 1, 0),
934
N( 0, 1),
935
N(-1, 0),
936
N( 0, -1)
937
in 2-d lattice N
938
sage: P1xP1 = ToricVariety(fan, "x s y t")
939
sage: P1xP1.inject_variables()
940
Defining x, s, y, t
941
sage: P1xP1._validate([x - y, x*s + y*t])
942
[x - y, x*s + y*t]
943
sage: P1xP1._validate([x + s])
944
Traceback (most recent call last):
945
...
946
ValueError: x + s is not homogeneous on
947
2-d toric variety covered by 4 affine patches!
948
"""
949
for p in polynomials:
950
if not self.is_homogeneous(p):
951
raise ValueError("%s is not homogeneous on %s!" % (p, self))
952
return polynomials
953
954
def affine_patch(self, i):
955
r"""
956
Return the ``i``-th affine patch of ``self``.
957
958
INPUT:
959
960
- ``i`` -- integer, index of a generating cone of the fan of ``self``.
961
962
OUTPUT:
963
964
- affine :class:`toric variety <ToricVariety_field>` corresponding to
965
the ``i``-th generating cone of the fan of ``self``.
966
967
The result is cached, so the ``i``-th patch is always the same object
968
in memory.
969
970
See also :meth:`affine_algebraic_patch`, which expresses the
971
patches as subvarieties of affine space instead.
972
973
EXAMPLES::
974
975
sage: fan = FaceFan(lattice_polytope.octahedron(2))
976
sage: P1xP1 = ToricVariety(fan, "x s y t")
977
sage: patch0 = P1xP1.affine_patch(0)
978
sage: patch0
979
2-d affine toric variety
980
sage: patch0.embedding_morphism()
981
Scheme morphism:
982
From: 2-d affine toric variety
983
To: 2-d toric variety covered by 4 affine patches
984
Defn: Defined on coordinates by sending [x : t] to
985
[x : 1 : 1 : t]
986
sage: patch1 = P1xP1.affine_patch(1)
987
sage: patch1.embedding_morphism()
988
Scheme morphism:
989
From: 2-d affine toric variety
990
To: 2-d toric variety covered by 4 affine patches
991
Defn: Defined on coordinates by sending [y : t] to
992
[1 : 1 : y : t]
993
sage: patch1 is P1xP1.affine_patch(1)
994
True
995
"""
996
i = int(i) # implicit type checking
997
try:
998
return self._affine_patches[i]
999
except AttributeError:
1000
self._affine_patches = dict()
1001
except KeyError:
1002
pass
1003
cone = self.fan().generating_cone(i)
1004
names = self.variable_names()
1005
# Number of "honest fan coordinates"
1006
n = self.fan().nrays()
1007
# Number of "torus factor coordinates"
1008
t = self._torus_factor_dim
1009
names = ([names[ray] for ray in cone.ambient_ray_indices()]
1010
+ list(names[n:]))
1011
patch = AffineToricVariety(cone, names, base_field=self.base_ring())
1012
embedding_coordinates = [1] * n
1013
for k, ray in enumerate(cone.ambient_ray_indices()):
1014
embedding_coordinates[ray] = patch.gen(k)
1015
if t > 0: # Passing "-0" gives unintended result
1016
embedding_coordinates.extend(patch.gens()[-t:])
1017
patch._embedding_morphism = patch.hom(embedding_coordinates, self)
1018
self._affine_patches[i] = patch
1019
return patch
1020
1021
def change_ring(self, F):
1022
r"""
1023
Return a toric variety over ``F`` and otherwise the same as ``self``.
1024
1025
INPUT:
1026
1027
- ``F`` -- field.
1028
1029
OUTPUT:
1030
1031
- :class:`toric variety <ToricVariety_field>` over ``F``.
1032
1033
.. NOTE::
1034
1035
There is no need to have any relation between ``F`` and the base
1036
field of ``self``. If you do want to have such a relation, use
1037
:meth:`base_extend` instead.
1038
1039
EXAMPLES::
1040
1041
sage: fan = FaceFan(lattice_polytope.octahedron(2))
1042
sage: P1xP1 = ToricVariety(fan)
1043
sage: P1xP1.base_ring()
1044
Rational Field
1045
sage: P1xP1_RR = P1xP1.change_ring(RR)
1046
sage: P1xP1_RR.base_ring()
1047
Real Field with 53 bits of precision
1048
sage: P1xP1_QQ = P1xP1_RR.change_ring(QQ)
1049
sage: P1xP1_QQ.base_ring()
1050
Rational Field
1051
sage: P1xP1_RR.base_extend(QQ)
1052
Traceback (most recent call last):
1053
...
1054
ValueError: no natural map from the base ring
1055
(=Real Field with 53 bits of precision)
1056
to R (=Rational Field)!
1057
"""
1058
if self.base_ring() == F:
1059
return self
1060
else:
1061
return ToricVariety(self.fan(), self.variable_names(),
1062
base_field=F)
1063
1064
def coordinate_ring(self):
1065
r"""
1066
Return the coordinate ring of ``self``.
1067
1068
For toric varieties this is the homogeneous coordinate ring (a.k.a.
1069
Cox's ring and total ring).
1070
1071
OUTPUT:
1072
1073
- polynomial ring.
1074
1075
EXAMPLES::
1076
1077
sage: fan = FaceFan(lattice_polytope.octahedron(2))
1078
sage: P1xP1 = ToricVariety(fan)
1079
sage: P1xP1.coordinate_ring()
1080
Multivariate Polynomial Ring in z0, z1, z2, z3
1081
over Rational Field
1082
1083
TESTS::
1084
1085
sage: R = toric_varieties.A1().coordinate_ring(); R
1086
Multivariate Polynomial Ring in z over Rational Field
1087
sage: type(R)
1088
<type 'sage.rings.polynomial.multi_polynomial_libsingular.MPolynomialRing_libsingular'>
1089
"""
1090
if "_coordinate_ring" not in self.__dict__:
1091
names = self.variable_names()
1092
self._coordinate_ring = PolynomialRing(self.base_ring(), len(names), names)
1093
return self._coordinate_ring
1094
1095
def embedding_morphism(self):
1096
r"""
1097
Return the default embedding morphism of ``self``.
1098
1099
Such a morphism is always defined for an affine patch of a toric
1100
variety (which is also a toric varieties itself).
1101
1102
OUTPUT:
1103
1104
- :class:`scheme morphism
1105
<sage.schemes.generic.morphism.SchemeMorphism_polynomial_toric_variety>`
1106
if the default embedding morphism was defined for ``self``,
1107
otherwise a ``ValueError`` exception is raised.
1108
1109
EXAMPLES::
1110
1111
sage: fan = FaceFan(lattice_polytope.octahedron(2))
1112
sage: P1xP1 = ToricVariety(fan, "x s y t")
1113
sage: P1xP1.embedding_morphism()
1114
Traceback (most recent call last):
1115
...
1116
ValueError: no default embedding was
1117
defined for this toric variety!
1118
sage: patch = P1xP1.affine_patch(0)
1119
sage: patch
1120
2-d affine toric variety
1121
sage: patch.embedding_morphism()
1122
Scheme morphism:
1123
From: 2-d affine toric variety
1124
To: 2-d toric variety covered by 4 affine patches
1125
Defn: Defined on coordinates by sending [x : t] to
1126
[x : 1 : 1 : t]
1127
"""
1128
try:
1129
return self._embedding_morphism
1130
except AttributeError:
1131
raise ValueError("no default embedding was defined for this "
1132
"toric variety!")
1133
1134
def fan(self, dim=None, codim=None):
1135
r"""
1136
Return the underlying fan of ``self`` or its cones.
1137
1138
INPUT:
1139
1140
- ``dim`` -- dimension of the requested cones;
1141
1142
- ``codim`` -- codimension of the requested cones.
1143
1144
OUTPUT:
1145
1146
- :class:`rational polyhedral fan
1147
<sage.geometry.fan.RationalPolyhedralFan>` if no parameters were
1148
given, :class:`tuple` of :class:`cones
1149
<sage.geometry.cone.ConvexRationalPolyhedralCone>` otherwise.
1150
1151
EXAMPLES::
1152
1153
sage: fan = FaceFan(lattice_polytope.octahedron(2))
1154
sage: P1xP1 = ToricVariety(fan)
1155
sage: P1xP1.fan()
1156
Rational polyhedral fan in 2-d lattice N
1157
sage: P1xP1.fan() is fan
1158
True
1159
sage: P1xP1.fan(1)[0]
1160
1-d cone of Rational polyhedral fan in 2-d lattice N
1161
"""
1162
return self._fan(dim, codim)
1163
1164
def inject_coefficients(self, scope=None, verbose=True):
1165
r"""
1166
Inject generators of the base field of ``self`` into ``scope``.
1167
1168
This function is useful if the base field is the field of rational
1169
functions.
1170
1171
INPUT:
1172
1173
- ``scope`` -- namespace (default: global, not just the scope from
1174
which this function was called);
1175
1176
- ``verbose`` -- if ``True`` (default), names of injected generators
1177
will be printed.
1178
1179
OUTPUT:
1180
1181
- none.
1182
1183
EXAMPLES::
1184
1185
sage: fan = FaceFan(lattice_polytope.octahedron(2))
1186
sage: F = QQ["a, b"].fraction_field()
1187
sage: P1xP1 = ToricVariety(fan, base_field=F)
1188
sage: P1xP1.inject_coefficients()
1189
Defining a, b
1190
1191
We check that we can use names ``a`` and ``b``, Trac #10498 is fixed::
1192
1193
sage: a + b
1194
a + b
1195
sage: a + b in P1xP1.coordinate_ring()
1196
True
1197
"""
1198
if scope is None:
1199
# scope = globals() does not work well, the above doctest fails.
1200
# Instead we "borrow" this code from sage.misc.misc.inject_variable
1201
depth = 0
1202
while True:
1203
scope = sys._getframe(depth).f_globals
1204
if (scope["__name__"] == "__main__"
1205
and scope["__package__"] is None):
1206
break
1207
depth += 1
1208
try:
1209
self.base_ring().inject_variables(scope, verbose)
1210
except AttributeError:
1211
pass
1212
1213
@cached_method
1214
def dimension_singularities(self):
1215
r"""
1216
Return the dimension of the singular set.
1217
1218
OUTPUT:
1219
1220
Integer. The dimension of the singular set of the toric
1221
variety. Often the singular set is a reducible subvariety, and
1222
this method will return the dimension of the
1223
largest-dimensional component.
1224
1225
Returns -1 if the toric variety is smooth.
1226
1227
EXAMPLES::
1228
1229
sage: toric_varieties.P4_11169().dimension_singularities()
1230
1
1231
sage: toric_varieties.Conifold().dimension_singularities()
1232
0
1233
sage: toric_varieties.P2().dimension_singularities()
1234
-1
1235
"""
1236
for codim in range(0,self.dimension()+1):
1237
if any(not cone.is_smooth() for cone in self.fan(codim)):
1238
return self.dimension() - codim
1239
return -1
1240
1241
def is_homogeneous(self, polynomial):
1242
r"""
1243
Check if ``polynomial`` is homogeneous.
1244
1245
The coordinate ring of a toric variety is multigraded by relations
1246
between generating rays of the underlying fan.
1247
1248
INPUT:
1249
1250
- ``polynomial`` -- polynomial in the coordinate ring of ``self`` or
1251
its quotient.
1252
1253
OUTPUT:
1254
1255
- ``True`` if ``polynomial`` is homogeneous and ``False`` otherwise.
1256
1257
EXAMPLES:
1258
1259
We will use the product of two projective lines with coordinates
1260
`(x, y)` for one and `(s, t)` for the other::
1261
1262
sage: fan = FaceFan(lattice_polytope.octahedron(2))
1263
sage: fan.rays()
1264
N( 1, 0),
1265
N( 0, 1),
1266
N(-1, 0),
1267
N( 0, -1)
1268
in 2-d lattice N
1269
sage: P1xP1 = ToricVariety(fan, "x s y t")
1270
sage: P1xP1.inject_variables()
1271
Defining x, s, y, t
1272
sage: P1xP1.is_homogeneous(x - y)
1273
True
1274
sage: P1xP1.is_homogeneous(x*s + y*t)
1275
True
1276
sage: P1xP1.is_homogeneous(x - t)
1277
False
1278
sage: P1xP1.is_homogeneous(1)
1279
True
1280
1281
Note that by homogeneous, we mean well-defined with respect to
1282
the homogeneous rescalings of ``self``. So a polynomial that you would
1283
usually not call homogeneous can be homogeneous if there are
1284
no homogeneous rescalings, for example::
1285
1286
sage: A1.<z> = toric_varieties.A1()
1287
sage: A1.is_homogeneous(z^3+z^7)
1288
True
1289
1290
Finally, the degree group is really the Chow group
1291
`A_{d-1}(X)` and can contain torsion. For example, take
1292
`\CC^2/\ZZ_2`. Here, the Chow group is `A_{d-1}(\CC^2/\ZZ_2) =
1293
\ZZ_2` and distinguishes even-degree homogeneous polynomials
1294
from odd-degree homogeneous polynomials::
1295
1296
sage: A2_Z2.<x,y> = toric_varieties.A2_Z2()
1297
sage: A2_Z2.is_homogeneous(x+y+x^3+y^5+x^3*y^4)
1298
True
1299
sage: A2_Z2.is_homogeneous(x^2+x*y+y^4+(x*y)^5+x^4*y^4)
1300
True
1301
sage: A2_Z2.is_homogeneous(x+y^2)
1302
False
1303
"""
1304
if '_homogeneous_degrees_group' not in self.__dict__:
1305
fan = self.fan()
1306
from sage.modules.free_module import FreeModule
1307
degrees_group = FreeModule(ZZ, fan.nrays()).quotient(
1308
fan.rays().matrix().columns())
1309
self._homogeneous_degrees_group = degrees_group
1310
degrees_group = self._homogeneous_degrees_group
1311
S = self.coordinate_ring()
1312
try:
1313
polynomial = S(polynomial)
1314
except TypeError:
1315
# Then it should be in the quotient corresponding to a subscheme
1316
polynomial = S(polynomial.lift())
1317
monomials = polynomial.monomials()
1318
if not monomials:
1319
return True
1320
degree = degrees_group(vector(ZZ,monomials[0].degrees()))
1321
for monomial in monomials:
1322
if degrees_group(vector(ZZ,monomial.degrees())) != degree:
1323
return False
1324
return True
1325
1326
def is_isomorphic(self, another):
1327
r"""
1328
Check if ``self`` is isomorphic to ``another``.
1329
1330
INPUT:
1331
1332
- ``another`` - :class:`toric variety <ToricVariety_field>`.
1333
1334
OUTPUT:
1335
1336
- ``True`` if ``self`` and ``another`` are isomorphic,
1337
``False`` otherwise.
1338
1339
EXAMPLES::
1340
1341
sage: fan1 = FaceFan(lattice_polytope.octahedron(2))
1342
sage: TV1 = ToricVariety(fan1)
1343
sage: fan2 = NormalFan(lattice_polytope.octahedron(2))
1344
sage: TV2 = ToricVariety(fan2)
1345
1346
Only the most trivial case is implemented so far::
1347
1348
sage: TV1.is_isomorphic(TV1)
1349
True
1350
sage: TV1.is_isomorphic(TV2)
1351
Traceback (most recent call last):
1352
...
1353
NotImplementedError:
1354
isomorphism check is not yet implemented!
1355
"""
1356
if self is another:
1357
return True
1358
if not is_ToricVariety(another):
1359
raise TypeError(
1360
"only another toric variety can be checked for isomorphism! "
1361
"Got %s" % another)
1362
raise NotImplementedError("isomorphism check is not yet implemented!")
1363
1364
def is_affine(self):
1365
r"""
1366
Check if ``self`` is an affine toric variety.
1367
1368
An affine toric variety is a toric variety whose fan is the
1369
face lattice of a single cone. See also
1370
:func:`AffineToricVariety`.
1371
1372
OUTPUT:
1373
1374
Boolean.
1375
1376
EXAMPLES::
1377
1378
sage: toric_varieties.A2().is_affine()
1379
True
1380
sage: toric_varieties.P1xA1().is_affine()
1381
False
1382
"""
1383
return self.fan().ngenerating_cones() == 1
1384
1385
def is_complete(self):
1386
r"""
1387
Check if ``self`` is complete.
1388
1389
OUTPUT:
1390
1391
- ``True`` if ``self`` is complete and ``False`` otherwise.
1392
1393
EXAMPLES::
1394
1395
sage: fan = FaceFan(lattice_polytope.octahedron(2))
1396
sage: P1xP1 = ToricVariety(fan)
1397
sage: P1xP1.is_complete()
1398
True
1399
sage: P1xP1.affine_patch(0).is_complete()
1400
False
1401
"""
1402
return self.fan().is_complete()
1403
1404
def is_orbifold(self):
1405
r"""
1406
Check if ``self`` has only quotient singularities.
1407
1408
A toric variety with at most orbifold singularities (in this
1409
sense) is often called a simplicial toric variety. In this
1410
package, we generally try to avoid this term since it mixes up
1411
differential geometry and cone terminology.
1412
1413
OUTPUT:
1414
1415
- ``True`` if ``self`` has at most quotient singularities by
1416
finite groups, ``False`` otherwise.
1417
1418
EXAMPLES::
1419
1420
sage: fan1 = FaceFan(lattice_polytope.octahedron(2))
1421
sage: P1xP1 = ToricVariety(fan1)
1422
sage: P1xP1.is_orbifold()
1423
True
1424
sage: fan2 = NormalFan(lattice_polytope.octahedron(3))
1425
sage: TV = ToricVariety(fan2)
1426
sage: TV.is_orbifold()
1427
False
1428
"""
1429
return self.fan().is_simplicial()
1430
1431
def is_smooth(self):
1432
r"""
1433
Check if ``self`` is smooth.
1434
1435
OUTPUT:
1436
1437
- ``True`` if ``self`` is smooth and ``False`` otherwise.
1438
1439
EXAMPLES::
1440
1441
sage: diamond = lattice_polytope.octahedron(2)
1442
sage: fan1 = FaceFan(diamond)
1443
sage: P1xP1 = ToricVariety(fan1)
1444
sage: P1xP1.is_smooth()
1445
True
1446
sage: fan2 = NormalFan(lattice_polytope.octahedron(2))
1447
sage: TV = ToricVariety(fan2)
1448
sage: TV.is_smooth()
1449
False
1450
"""
1451
return self.fan().is_smooth()
1452
1453
@cached_method
1454
def Kaehler_cone(self):
1455
r"""
1456
Return the closure of the Kähler cone of ``self``.
1457
1458
OUTPUT:
1459
1460
- :class:`cone <sage.geometry.cone.ConvexRationalPolyhedralCone>`.
1461
1462
.. NOTE::
1463
1464
This cone sits in the rational divisor class group of ``self`` and
1465
the choice of coordinates agrees with
1466
:meth:`rational_class_group`.
1467
1468
EXAMPLES::
1469
1470
sage: P1xP1 = toric_varieties.P1xP1()
1471
sage: Kc = P1xP1.Kaehler_cone()
1472
sage: Kc
1473
2-d cone in 2-d lattice
1474
sage: Kc.rays()
1475
Divisor class [0, 1],
1476
Divisor class [1, 0]
1477
in Basis lattice of The toric rational divisor class group
1478
of a 2-d CPR-Fano toric variety covered by 4 affine patches
1479
sage: [ divisor_class.lift() for divisor_class in Kc.rays() ]
1480
[V(x), V(s)]
1481
sage: Kc.lattice()
1482
Basis lattice of The toric rational divisor class group of a
1483
2-d CPR-Fano toric variety covered by 4 affine patches
1484
"""
1485
fan = self.fan()
1486
GT = fan.Gale_transform().columns()
1487
from sage.schemes.toric.divisor import \
1488
ToricRationalDivisorClassGroup_basis_lattice
1489
L = ToricRationalDivisorClassGroup_basis_lattice(
1490
self.rational_class_group())
1491
n = fan.nrays()
1492
K = None
1493
for cone in fan:
1494
sigma = Cone([GT[i] for i in range(n)
1495
if i not in cone.ambient_ray_indices()],
1496
lattice=L)
1497
K = K.intersection(sigma) if K is not None else sigma
1498
return K
1499
1500
@cached_method
1501
def Mori_cone(self):
1502
r"""
1503
Returns the Mori cone of ``self``.
1504
1505
OUTPUT:
1506
1507
- :class:`cone <sage.geometry.cone.ConvexRationalPolyhedralCone>`.
1508
1509
.. NOTE::
1510
1511
* The Mori cone is dual to the Kähler cone.
1512
1513
* We think of the Mori cone as living inside the row span of the
1514
Gale transform matrix (computed by
1515
``self.fan().Gale_transform()``).
1516
1517
* The points in the Mori cone are the effective curves in the
1518
variety.
1519
1520
* The ``i``-th entry in each Mori vector is the intersection
1521
number of the curve corresponding to the generator of the
1522
``i``-th ray of the fan with the corresponding divisor class.
1523
The very last entry is associated to the orgin of the fan
1524
lattice.
1525
1526
* The Mori vectors are also known as the gauged linear sigma model
1527
charge vectors.
1528
1529
EXAMPLES::
1530
1531
sage: P4_11169 = toric_varieties.P4_11169_resolved()
1532
sage: P4_11169.Mori_cone()
1533
2-d cone in 7-d lattice
1534
sage: P4_11169.Mori_cone().rays()
1535
(0, 0, 1, 1, 1, -3, 0),
1536
(3, 2, 0, 0, 0, 1, -6)
1537
in Ambient free module of rank 7
1538
over the principal ideal domain Integer Ring
1539
"""
1540
# Ideally, self.Kaehler_cone().dual() should be it, but
1541
# so far this is not the case.
1542
rays = (ray * self._fan.Gale_transform()
1543
for ray in self.Kaehler_cone().dual().rays())
1544
return Cone(rays, lattice=ZZ**(self._fan.nrays()+1))
1545
1546
def plot(self, **options):
1547
r"""
1548
Plot ``self``, i.e. the corresponding fan.
1549
1550
INPUT:
1551
1552
- any options for toric plots (see :func:`toric_plotter.options
1553
<sage.geometry.toric_plotter.options>`), none are mandatory.
1554
1555
OUTPUT:
1556
1557
- a plot.
1558
1559
.. NOTE::
1560
1561
The difference between ``X.plot()`` and ``X.fan().plot()`` is that
1562
in the first case default ray labels correspond to variables of
1563
``X``.
1564
1565
EXAMPLES::
1566
1567
sage: X = toric_varieties.Cube_deformation(4)
1568
sage: X.plot()
1569
"""
1570
if "ray_label" not in options:
1571
gens = self.coordinate_ring().gens()
1572
if self.fan().lattice().degree() <= 2:
1573
options["ray_label"] = ["$%s$" % latex(z) for z in gens]
1574
else:
1575
options["ray_label"] = [str(z) for z in gens]
1576
return self.fan().plot(**options)
1577
1578
def rational_class_group(self):
1579
r"""
1580
Return the rational divisor class group of ``self``.
1581
1582
Let `X` be a toric variety.
1583
1584
The **Weil divisor class group** `\mathop{Cl}(X)` is a finitely
1585
generated abelian group and can contain torsion. Its rank equals the
1586
number of rays in the fan of `X` minus the dimension of `X`.
1587
1588
The **rational divisor class group** is
1589
`\mathop{Cl}(X) \otimes_\ZZ \QQ` and never includes torsion. If `X` is
1590
*smooth*, this equals the **Picard group** of `X`, whose elements are
1591
the isomorphism classes of line bundles on `X`. The group law (which
1592
we write as addition) is the tensor product of the line bundles. The
1593
Picard group of a toric variety is always torsion-free.
1594
1595
OUTPUT:
1596
1597
- :class:`rational divisor class group
1598
<sage.schemes.toric.divisor.ToricRationalDivisorClassGroup>`.
1599
1600
.. NOTE::
1601
1602
* Coordinates correspond to the rows of
1603
``self.fan().gale_transform()``.
1604
1605
* :meth:`Kaehler_cone` yields a cone in this group.
1606
1607
EXAMPLES::
1608
1609
sage: fan = FaceFan(lattice_polytope.octahedron(2))
1610
sage: P1xP1 = ToricVariety(fan)
1611
sage: P1xP1.rational_class_group()
1612
The toric rational divisor class group
1613
of a 2-d toric variety covered by 4 affine patches
1614
"""
1615
from sage.schemes.toric.divisor import ToricRationalDivisorClassGroup
1616
return ToricRationalDivisorClassGroup(self)
1617
1618
def Chow_group(self, base_ring=ZZ):
1619
r"""
1620
Return the toric Chow group.
1621
1622
INPUT:
1623
1624
- ``base_ring`` -- either ``ZZ`` (default) or ``QQ``. The
1625
coefficient ring of the Chow group.
1626
1627
OUTPUT:
1628
1629
A :class:`sage.schemes.toric.chow_group.ChowGroup_class`
1630
1631
EXAMPLES::
1632
1633
sage: A = toric_varieties.P2().Chow_group(); A
1634
Chow group of 2-d CPR-Fano toric variety covered by 3 affine patches
1635
sage: A.gens()
1636
(( 1 | 0 | 0 ), ( 0 | 1 | 0 ), ( 0 | 0 | 1 ))
1637
"""
1638
from sage.schemes.toric.chow_group import ChowGroup
1639
return ChowGroup(self,base_ring)
1640
1641
def cartesian_product(self, other,
1642
coordinate_names=None, coordinate_indices=None):
1643
r"""
1644
Return the Cartesian product of ``self`` with ``other``.
1645
1646
INPUT:
1647
1648
- ``other`` -- a :class:`toric variety <ToricVariety_field>`;
1649
1650
- ``coordinate_names`` -- names of variables for the coordinate ring,
1651
see :func:`normalize_names` for acceptable formats. If not given,
1652
indexed variable names will be created automatically;
1653
1654
- ``coordinate_indices`` -- list of integers, indices for indexed
1655
variables. If not given, the index of each variable will coincide
1656
with the index of the corresponding ray of the fan.
1657
1658
OUTPUT:
1659
1660
-- a :class:`toric variety <ToricVariety_field>`.
1661
1662
EXAMPLES::
1663
1664
sage: P1 = ToricVariety(Fan([Cone([(1,)]), Cone([(-1,)])]))
1665
sage: P1xP1 = P1.cartesian_product(P1); P1xP1
1666
2-d toric variety covered by 4 affine patches
1667
sage: P1xP1.fan().rays()
1668
N+N(-1, 0),
1669
N+N( 1, 0),
1670
N+N( 0, -1),
1671
N+N( 0, 1)
1672
in 2-d lattice N+N
1673
"""
1674
return ToricVariety(self.fan().cartesian_product(other.fan()),
1675
coordinate_names, coordinate_indices,
1676
base_field=self.base_ring())
1677
1678
def resolve(self, **kwds):
1679
r"""
1680
Construct a toric variety whose fan subdivides the fan of ``self``.
1681
1682
The name of this function reflects the fact that usually such
1683
subdivisions are done for resolving singularities of the original
1684
variety.
1685
1686
INPUT:
1687
1688
This function accepts only keyword arguments, none of which are
1689
mandatory.
1690
1691
- ``coordinate_names`` -- names for coordinates of the new variety. If
1692
not given, will be constructed from the coordinate names of ``self``
1693
and necessary indexed ones. See :func:`normalize_names` for the
1694
description of acceptable formats;
1695
1696
- ``coordinate_indices`` -- coordinate indices which should be used
1697
for indexed variables of the new variety;
1698
1699
- all other arguments will be passed to
1700
:meth:`~sage.geometry.fan.RationalPolyhedralFan.subdivide` method of
1701
the underlying :class:`rational polyhedral fan
1702
<sage.geometry.fan.RationalPolyhedralFan>`, see its documentation
1703
for the available options.
1704
1705
OUTPUT:
1706
1707
- :class:`toric variety <ToricVariety_field>`.
1708
1709
EXAMPLES:
1710
1711
First we will "manually" resolve a simple orbifold singularity::
1712
1713
sage: cone = Cone([(1,1), (-1,1)])
1714
sage: fan = Fan([cone])
1715
sage: TV = ToricVariety(fan)
1716
sage: TV.is_smooth()
1717
False
1718
sage: TV_res = TV.resolve(new_rays=[(0,1)])
1719
sage: TV_res.is_smooth()
1720
True
1721
sage: TV_res.fan().rays()
1722
N( 1, 1),
1723
N(-1, 1),
1724
N( 0, 1)
1725
in 2-d lattice N
1726
sage: [cone.ambient_ray_indices() for cone in TV_res.fan()]
1727
[(0, 2), (1, 2)]
1728
1729
Now let's "automatically" partially resolve a more complicated fan::
1730
1731
sage: fan = NormalFan(lattice_polytope.octahedron(3))
1732
sage: TV = ToricVariety(fan)
1733
sage: TV.is_smooth()
1734
False
1735
sage: TV.is_orbifold()
1736
False
1737
sage: TV.fan().nrays()
1738
8
1739
sage: TV.fan().ngenerating_cones()
1740
6
1741
sage: TV_res = TV.resolve(make_simplicial=True)
1742
sage: TV_res.is_smooth()
1743
False
1744
sage: TV_res.is_orbifold()
1745
True
1746
sage: TV_res.fan().nrays()
1747
8
1748
sage: TV_res.fan().ngenerating_cones()
1749
12
1750
sage: TV.gens()
1751
(z0, z1, z2, z3, z4, z5, z6, z7)
1752
sage: TV_res.gens()
1753
(z0, z1, z2, z3, z4, z5, z6, z7)
1754
sage: TV_res = TV.resolve(coordinate_names="x+",
1755
... make_simplicial=True)
1756
sage: TV_res.gens()
1757
(x0, x1, x2, x3, x4, x5, x6, x7)
1758
"""
1759
# If you are changing this function, check out resolve in Fano toric
1760
# varieties to see if it should be changed too
1761
#
1762
# Currently the resolution of fans works for full-dimensional ones
1763
# only, so there is no point to deal with the general case here, since
1764
# we will not be able to check that it works.
1765
coordinate_names = kwds.pop("coordinate_names", None)
1766
coordinate_indices = kwds.pop("coordinate_indices", None)
1767
fan = self.fan()
1768
if fan.dim() != fan.lattice_dim():
1769
raise NotImplementedError("resolution of toric varieties with "
1770
"torus factors is not yet implemented!")
1771
# When it is implemented, should be careful with the torus factor
1772
rfan = fan.subdivide(**kwds)
1773
if coordinate_names is None:
1774
coordinate_names = list(self.variable_names())
1775
if coordinate_indices is None:
1776
coordinate_indices = range(fan.nrays(), rfan.nrays())
1777
else:
1778
coordinate_indices = coordinate_indices[fan.nrays():]
1779
coordinate_names.extend(normalize_names(
1780
ngens=rfan.nrays() - fan.nrays(),
1781
indices=coordinate_indices,
1782
prefix=self._coordinate_prefix))
1783
coordinate_names.append(self._coordinate_prefix + "+")
1784
resolution = ToricVariety(rfan, coordinate_names=coordinate_names,
1785
coordinate_indices=coordinate_indices,
1786
base_field=self.base_ring())
1787
R = self.coordinate_ring()
1788
R_res = resolution.coordinate_ring()
1789
resolution_map = resolution.hom(R.hom(R_res.gens()[:R.ngens()]), self)
1790
resolution._resolution_map = resolution_map
1791
# The above map does not have (yet) public methods to access it.
1792
# While this map is defined correctly, base classes of schemes and
1793
# morphisms do not treat it as they should. The plan is to fix this
1794
# situation soon and to be able to use this map!
1795
return resolution
1796
1797
def resolve_to_orbifold(self, **kwds):
1798
r"""
1799
Construct an orbifold whose fan subdivides the fan of ``self``.
1800
1801
It is a synonym for :meth:`resolve` with ``make_simplicial=True``
1802
option.
1803
1804
INPUT:
1805
1806
- this function accepts only keyword arguments. See :meth:`resolve`
1807
for documentation.
1808
1809
OUTPUT:
1810
1811
- :class:`toric variety <ToricVariety_field>`.
1812
1813
EXAMPLES::
1814
1815
sage: fan = NormalFan(lattice_polytope.octahedron(3))
1816
sage: TV = ToricVariety(fan)
1817
sage: TV.is_orbifold()
1818
False
1819
sage: TV.fan().nrays()
1820
8
1821
sage: TV.fan().ngenerating_cones()
1822
6
1823
sage: TV_res = TV.resolve_to_orbifold()
1824
sage: TV_res.is_orbifold()
1825
True
1826
sage: TV_res.fan().nrays()
1827
8
1828
sage: TV_res.fan().ngenerating_cones()
1829
12
1830
"""
1831
return self.resolve(make_simplicial=True, **kwds)
1832
1833
def subscheme(self, polynomials):
1834
r"""
1835
Return the subscheme of ``self`` defined by ``polynomials``.
1836
1837
INPUT:
1838
1839
- ``polynomials`` -- list of polynomials in the coordinate ring of
1840
``self``.
1841
1842
OUTPUT:
1843
1844
- :class:`subscheme of a toric variety
1845
<sage.schemes.generic.algebraic_scheme.AlgebraicScheme_subscheme_toric>`.
1846
1847
EXAMPLES:
1848
1849
We will construct a subscheme of the product of two projective lines
1850
with coordinates `(x, y)` for one and `(s, t)` for the other::
1851
1852
sage: fan = FaceFan(lattice_polytope.octahedron(2))
1853
sage: fan.rays()
1854
N( 1, 0),
1855
N( 0, 1),
1856
N(-1, 0),
1857
N( 0, -1)
1858
in 2-d lattice N
1859
sage: P1xP1 = ToricVariety(fan, "x s y t")
1860
sage: P1xP1.inject_variables()
1861
Defining x, s, y, t
1862
sage: X = P1xP1.subscheme([x*s + y*t, x^3+y^3])
1863
sage: X
1864
Closed subscheme of 2-d toric variety
1865
covered by 4 affine patches defined by:
1866
x*s + y*t,
1867
x^3 + y^3
1868
sage: X.defining_polynomials()
1869
(x*s + y*t, x^3 + y^3)
1870
sage: X.defining_ideal()
1871
Ideal (x*s + y*t, x^3 + y^3)
1872
of Multivariate Polynomial Ring in x, s, y, t
1873
over Rational Field
1874
sage: X.base_ring()
1875
Rational Field
1876
sage: X.base_scheme()
1877
Spectrum of Rational Field
1878
sage: X.structure_morphism()
1879
Scheme morphism:
1880
From: Closed subscheme of
1881
2-d toric variety covered by 4 affine patches defined by:
1882
x*s + y*t,
1883
x^3 + y^3
1884
To: Spectrum of Rational Field
1885
Defn: Structure map
1886
"""
1887
from sage.schemes.generic.algebraic_scheme import \
1888
AlgebraicScheme_subscheme_toric, AlgebraicScheme_subscheme_affine_toric
1889
if self.is_affine():
1890
return AlgebraicScheme_subscheme_affine_toric(self, polynomials)
1891
else:
1892
return AlgebraicScheme_subscheme_toric(self, polynomials)
1893
1894
def Stanley_Reisner_ideal(self):
1895
"""
1896
Return the Stanley-Reisner ideal.
1897
1898
OUTPUT:
1899
1900
- The Stanley-Reisner ideal in the polynomial ring over
1901
`\QQ` generated by the homogeneous coordinates.
1902
1903
EXAMPLES::
1904
1905
sage: fan = Fan([[0,1,3],[3,4],[2,0],[1,2,4]], [(-3, -2, 1), (0, 0, 1), (3, -2, 1), (-1, -1, 1), (1, -1, 1)])
1906
sage: X = ToricVariety(fan, coordinate_names='A B C D E', base_field=GF(5))
1907
sage: SR = X.Stanley_Reisner_ideal(); SR
1908
Ideal (A*E, C*D, A*B*C, B*D*E) of Multivariate Polynomial Ring in A, B, C, D, E over Rational Field
1909
"""
1910
if "_SR" not in self.__dict__:
1911
R = PolynomialRing(QQ, self.variable_names())
1912
self._SR = self._fan.Stanley_Reisner_ideal(R)
1913
return self._SR
1914
1915
def linear_equivalence_ideal(self):
1916
"""
1917
Return the ideal generated by linear relations
1918
1919
OUTPUT:
1920
1921
- The ideal generated by the linear relations of the rays in
1922
the polynomial ring over `\QQ` generated by the homogeneous
1923
coordinates.
1924
1925
EXAMPLES::
1926
1927
sage: fan = Fan([[0,1,3],[3,4],[2,0],[1,2,4]], [(-3, -2, 1), (0, 0, 1), (3, -2, 1), (-1, -1, 1), (1, -1, 1)])
1928
sage: X = ToricVariety(fan, coordinate_names='A B C D E', base_field=GF(5))
1929
sage: lin = X.linear_equivalence_ideal(); lin
1930
Ideal (-3*A + 3*C - D + E, -2*A - 2*C - D - E, A + B + C + D + E) of Multivariate Polynomial Ring in A, B, C, D, E over Rational Field
1931
"""
1932
if "_linear_equivalence_ideal" not in self.__dict__:
1933
R = PolynomialRing(QQ, self.variable_names())
1934
self._linear_equivalence_ideal = self._fan.linear_equivalence_ideal(R)
1935
return self._linear_equivalence_ideal
1936
1937
@cached_method
1938
def cohomology_ring(self):
1939
r"""
1940
Return the cohomology ring of the toric variety.
1941
1942
OUTPUT:
1943
1944
- If the toric variety is is over `\CC` and has at most finite
1945
orbifold singularities: `H^\bullet(X,\QQ)` as a polynomial
1946
quotient ring.
1947
1948
- Other cases are not handled yet.
1949
1950
.. NOTE::
1951
1952
- Toric varieties over any field of characteristic 0 are
1953
treated as if they were varieties over `\CC`.
1954
1955
- The integral cohomology of smooth toric varieties is
1956
torsion-free, so in this case there is no loss of
1957
information when going to rational coefficients.
1958
1959
- ``self.cohomology_ring().gen(i)`` is the divisor class corresponding to
1960
the ``i``-th ray of the fan.
1961
1962
EXAMPLES::
1963
1964
sage: X = toric_varieties.dP6()
1965
sage: X.cohomology_ring()
1966
Rational cohomology ring of a 2-d CPR-Fano toric variety covered by 6 affine patches
1967
sage: X.cohomology_ring().defining_ideal()
1968
Ideal (-u - y + z + w, x - y - v + w, x*y, x*v, x*z, u*v, u*z, u*w, y*z, y*w, v*w) of Multivariate Polynomial Ring in x, u, y, v, z, w over Rational Field
1969
sage: X.cohomology_ring().defining_ideal().ring()
1970
Multivariate Polynomial Ring in x, u, y, v, z, w over Rational Field
1971
sage: X.variable_names()
1972
('x', 'u', 'y', 'v', 'z', 'w')
1973
sage: X.cohomology_ring().gens()
1974
([y + v - w], [-y + z + w], [y], [v], [z], [w])
1975
1976
TESTS:
1977
1978
The cohomology ring is a circular reference that is
1979
potentially troublesome on unpickling, see :trac:`15050`
1980
and :trac:`15149` ::
1981
1982
sage: variety = toric_varieties.P(1)
1983
sage: a = [variety.cohomology_ring(), variety.cohomology_basis(), variety.volume_class()]
1984
sage: b = [variety.Todd_class(), variety.Chern_class(), variety.Chern_character(), variety.Kaehler_cone(), variety.Mori_cone()]
1985
sage: loads(dumps(variety)) == variety
1986
True
1987
"""
1988
if self.base_ring().characteristic()>0:
1989
raise NotImplementedError('Only characteristic 0 base fields '
1990
'are implemented.')
1991
return CohomologyRing(self)
1992
1993
@cached_method
1994
def cohomology_basis(self, d=None):
1995
r"""
1996
Return a basis for the cohomology of the toric variety.
1997
1998
INPUT:
1999
2000
- ``d`` (optional) -- integer.
2001
2002
OUTPUT:
2003
2004
- Without the optional argument, a list whose d-th entry is a
2005
basis for `H^{2d}(X,\QQ)`
2006
2007
- If the argument is an integer ``d``, returns basis for
2008
`H^{2d}(X,\QQ)`
2009
2010
EXAMPLES::
2011
2012
sage: X = toric_varieties.dP8()
2013
sage: X.cohomology_basis()
2014
(([1],), ([z], [y]), ([y*z],))
2015
sage: X.cohomology_basis(1)
2016
([z], [y])
2017
sage: X.cohomology_basis(dimension(X))[0] == X.volume_class()
2018
True
2019
"""
2020
if d!=None:
2021
return self.cohomology_basis()[d]
2022
2023
H = self.cohomology_ring()
2024
# Make an empty list for each d-piece
2025
basis = [[] for d in range(self.dimension() + 1)]
2026
# Distribute basis elements into d-pieces
2027
for x in H.defining_ideal().normal_basis():
2028
basis[x.total_degree()].append(x)
2029
# Convert list of lists of polynomials to
2030
# tuple of tuples of cohomology classes
2031
return tuple(tuple(H(x) for x in dbasis)
2032
for dbasis in basis)
2033
2034
@cached_method
2035
def volume_class(self):
2036
r"""
2037
Return the cohomology class of the volume form on the toric
2038
variety.
2039
2040
Note that we are using cohomology with compact supports. If
2041
the variety is non-compact this is dual to homology without
2042
any support condition. In particular, for non-compact
2043
varieties the volume form `\mathrm{dVol}=\wedge_i(dx_i \wedge
2044
dy_i)` does not define a (non-zero) cohomology class.
2045
2046
OUTPUT:
2047
2048
A :class:`CohomologyClass`. If it exists, it is the class of
2049
the (properly normalized) volume form, that is, it is the
2050
Poincare dual of a single point. If it does not exist, a
2051
``ValueError`` is raised.
2052
2053
EXAMPLES::
2054
2055
sage: P2 = toric_varieties.P2()
2056
sage: P2.volume_class()
2057
[z^2]
2058
2059
sage: A2_Z2 = toric_varieties.A2_Z2()
2060
sage: A2_Z2.volume_class()
2061
Traceback (most recent call last):
2062
...
2063
ValueError: Volume class does not exist.
2064
2065
If none of the maximal cones is smooth things get more
2066
tricky. In this case no torus-fixed point is smooth. If we
2067
want to count an ordinary point as `1`, then a `G`-orbifold
2068
point needs to count as `\frac{1}{|G|}`. For example, take
2069
`\mathbb{P}^1\times\mathbb{P}^1` with inhomogeneous
2070
coordinates `(t,y)`. Take the quotient by the action
2071
`(t,y)\mapsto (-t,-y)`. The `\ZZ_2`-invariant Weil divisors
2072
`\{t=0\}` and `\{y=0\}` intersect in a `\ZZ_2`-fixed point, so
2073
they ought to have intersection number `\frac{1}{2}`. This
2074
means that the cohomology class `[t] \cap [y]` should be
2075
`\frac{1}{2}` times the volume class. Note that this is
2076
different from the volume normalization chosen in
2077
[Schubert]_::
2078
2079
sage: P1xP1_Z2 = toric_varieties.P1xP1_Z2()
2080
sage: Dt = P1xP1_Z2.divisor(1); Dt
2081
V(t)
2082
sage: Dy = P1xP1_Z2.divisor(3); Dy
2083
V(y)
2084
sage: P1xP1_Z2.volume_class()
2085
[2*t*y]
2086
2087
sage: HH = P1xP1_Z2.cohomology_ring()
2088
sage: HH(Dt) * HH(Dy) == 1/2 * P1xP1_Z2.volume_class()
2089
True
2090
2091
The fractional coefficients are also necessary to match the
2092
normalization in the rational Chow group for simplicial toric
2093
varieties::
2094
2095
sage: A = P1xP1_Z2.Chow_group(QQ)
2096
sage: A(Dt).intersection_with_divisor(Dy).count_points()
2097
1/2
2098
2099
REFERENCES:
2100
2101
.. [Schubert]
2102
Sheldon Katz and Stein Arild Stromme,
2103
A Maple package for intersection theory and enumerative geometry.
2104
"""
2105
if not self.is_orbifold():
2106
raise NotImplementedError('Cohomology computations are only '
2107
'implemented for orbifolds.')
2108
HH = self.cohomology_ring()
2109
dim = self.dimension_relative()
2110
dVol = HH(self.fan().generating_cone(0)).part_of_degree(dim)
2111
if dVol.is_zero():
2112
raise ValueError, 'Volume class does not exist.'
2113
return dVol
2114
2115
def integrate(self, cohomology_class):
2116
"""
2117
Integrate a cohomology class over the toric variety.
2118
2119
INPUT:
2120
2121
- ``cohomology_class`` -- A cohomology class given as a
2122
polynomial in ``self.cohomology_ring()``
2123
2124
OUTPUT:
2125
2126
The integral of the cohomology class over the variety. The
2127
volume normalization is given by :meth:`volume_class`, that
2128
is, ``self.integrate(self.volume_class())`` is always one (if
2129
the volume class exists).
2130
2131
EXAMPLES::
2132
2133
sage: dP6 = toric_varieties.dP6()
2134
sage: HH = dP6.cohomology_ring()
2135
sage: D = [ HH(c) for c in dP6.fan(dim=1) ]
2136
sage: matrix([ [ D[i]*D[j] for i in range(0,6) ] for j in range(0,6) ])
2137
[ [w^2] [-w^2] [0] [0] [0] [-w^2]]
2138
[[-w^2] [w^2] [-w^2] [0] [0] [0]]
2139
[ [0] [-w^2] [w^2] [-w^2] [0] [0]]
2140
[ [0] [0] [-w^2] [w^2] [-w^2] [0]]
2141
[ [0] [0] [0] [-w^2] [w^2] [-w^2]]
2142
[[-w^2] [0] [0] [0] [-w^2] [w^2]]
2143
sage: matrix([ [ dP6.integrate(D[i]*D[j]) for i in range(0,6) ] for j in range(0,6) ])
2144
[-1 1 0 0 0 1]
2145
[ 1 -1 1 0 0 0]
2146
[ 0 1 -1 1 0 0]
2147
[ 0 0 1 -1 1 0]
2148
[ 0 0 0 1 -1 1]
2149
[ 1 0 0 0 1 -1]
2150
2151
If the toric variety is an orbifold, the intersection numbers
2152
are usually fractional::
2153
2154
sage: P2_123 = toric_varieties.P2_123()
2155
sage: HH = P2_123.cohomology_ring()
2156
sage: D = [ HH(c) for c in P2_123.fan(dim=1) ]
2157
sage: matrix([ [ P2_123.integrate(D[i]*D[j]) for i in range(0,3) ] for j in range(0,3) ])
2158
[2/3 1 1/3]
2159
[ 1 3/2 1/2]
2160
[1/3 1/2 1/6]
2161
sage: A = P2_123.Chow_group(QQ)
2162
sage: matrix([ [ A(P2_123.divisor(i))
2163
... .intersection_with_divisor(P2_123.divisor(j))
2164
... .count_points() for i in range(0,3) ] for j in range(0,3) ])
2165
[2/3 1 1/3]
2166
[ 1 3/2 1/2]
2167
[1/3 1/2 1/6]
2168
"""
2169
assert self.is_complete(), "Can only integrate over compact varieties."
2170
top_form = cohomology_class.part_of_degree(self.dimension())
2171
if top_form.is_zero(): return 0
2172
return top_form.lc() / self.volume_class().lc()
2173
2174
@cached_method
2175
def Chern_class(self, deg=None):
2176
"""
2177
Return Chern classes of the (tangent bundle of the) toric variety.
2178
2179
INPUT:
2180
2181
- ``deg`` -- integer (optional). The degree of the Chern class.
2182
2183
OUTPUT:
2184
2185
- If the degree is specified, the ``deg``-th Chern class.
2186
2187
- If no degree is specified, the total Chern class.
2188
2189
REFERENCES:
2190
2191
..
2192
2193
http://en.wikipedia.org/wiki/Chern_class
2194
2195
EXAMPLES::
2196
2197
sage: X = toric_varieties.dP6()
2198
sage: X.Chern_class()
2199
[-6*w^2 + y + 2*v + 2*z + w + 1]
2200
sage: X.c()
2201
[-6*w^2 + y + 2*v + 2*z + w + 1]
2202
sage: X.c(1)
2203
[y + 2*v + 2*z + w]
2204
sage: X.c(2)
2205
[-6*w^2]
2206
sage: X.integrate( X.c(2) )
2207
6
2208
sage: X.integrate( X.c(2) ) == X.Euler_number()
2209
True
2210
"""
2211
assert self.is_orbifold(), "Requires the toric variety to be an orbifold."
2212
c = prod([ 1+self.cohomology_ring().gen(i) for i in range(0,self._fan.nrays()) ])
2213
if deg==None:
2214
return c
2215
else:
2216
return c.part_of_degree(deg)
2217
2218
@cached_method
2219
def Chern_character(self, deg=None):
2220
"""
2221
Return the Chern character (of the tangent bundle) of the toric
2222
variety.
2223
2224
INPUT:
2225
2226
- ``deg`` -- integer (optional). The degree of the Chern
2227
character.
2228
2229
OUTPUT:
2230
2231
- If the degree is specified, the degree-``deg`` part of the
2232
Chern character.
2233
2234
- If no degree is specified, the total Chern character.
2235
2236
REFERENCES:
2237
2238
..
2239
2240
http://en.wikipedia.org/wiki/Chern_character#The_Chern_character
2241
2242
EXAMPLES::
2243
2244
sage: dP6 = toric_varieties.dP6()
2245
sage: dP6.Chern_character()
2246
[3*w^2 + y + 2*v + 2*z + w + 2]
2247
sage: dP6.ch()
2248
[3*w^2 + y + 2*v + 2*z + w + 2]
2249
sage: dP6.ch(1) == dP6.c(1)
2250
True
2251
"""
2252
assert self.is_orbifold(), "Requires the toric variety to be an orbifold."
2253
n_rels = self._fan.nrays() - self.dimension()
2254
ch = sum([ self.cohomology_ring().gen(i).exp()
2255
for i in range(0,self._fan.nrays()) ]) - n_rels
2256
if deg==None:
2257
return ch
2258
else:
2259
return ch.part_of_degree(deg)
2260
2261
@cached_method
2262
def Todd_class(self, deg=None):
2263
"""
2264
Return the Todd class (of the tangent bundle) of the toric variety.
2265
2266
INPUT:
2267
2268
- ``deg`` -- integer (optional). The desired degree part.
2269
2270
OUTPUT:
2271
2272
- If the degree is specified, the degree-``deg`` part of the
2273
Todd class.
2274
2275
- If no degree is specified, the total Todd class.
2276
2277
REFERENCES:
2278
2279
..
2280
2281
http://en.wikipedia.org/wiki/Todd_class
2282
2283
EXAMPLES::
2284
2285
sage: dP6 = toric_varieties.dP6()
2286
sage: dP6.Todd_class()
2287
[-w^2 + 1/2*y + v + z + 1/2*w + 1]
2288
sage: dP6.Td()
2289
[-w^2 + 1/2*y + v + z + 1/2*w + 1]
2290
sage: dP6.integrate( dP6.Td() )
2291
1
2292
"""
2293
Td = QQ(1)
2294
if self.dimension() >= 1:
2295
c1 = self.Chern_class(1)
2296
Td += QQ(1)/2 * c1
2297
if self.dimension() >= 2:
2298
c2 = self.Chern_class(2)
2299
Td += QQ(1)/12 * (c1**2 + c2)
2300
if self.dimension() >= 3:
2301
Td += QQ(1)/24 * c1*c2
2302
if self.dimension() >= 4:
2303
c3 = self.Chern_class(3)
2304
c4 = self.Chern_class(4)
2305
Td += -QQ(1)/720 * (c1**4 -4*c1**2*c2 -3*c2**2 -c1*c3 +c4)
2306
if self.dimension() >= 5:
2307
raise NotImplementedError('Todd class is currently only implemented up to degree 4')
2308
if deg==None:
2309
return Td
2310
else:
2311
return Td.part_of_degree(deg)
2312
2313
c = Chern_class
2314
ch = Chern_character
2315
Td = Todd_class
2316
2317
def Euler_number(self):
2318
"""
2319
Return the topological Euler number of the toric variety.
2320
2321
Sometimes, this is also called the Euler
2322
characteristic. :meth:`chi` is a synonym for
2323
:meth:`Euler_number`.
2324
2325
REFERENCES:
2326
2327
..
2328
2329
http://en.wikipedia.org/wiki/Euler_characteristic
2330
2331
EXAMPLES::
2332
2333
sage: P1xP1 = toric_varieties.P1xP1()
2334
sage: P1xP1.Euler_number()
2335
4
2336
sage: P1xP1.chi()
2337
4
2338
"""
2339
if "_chi" not in self.__dict__:
2340
if self.is_complete():
2341
chi = self.integrate(self.Chern_class())
2342
else:
2343
chi=0
2344
H = self.cohomology_basis()
2345
for d in range(0, self.dimension()+1):
2346
chi += (-1)**d * len(H[d])
2347
self._chi = chi
2348
return self._chi
2349
2350
chi = Euler_number
2351
2352
def K(self):
2353
r"""
2354
Returns the canonical divisor of the toric variety.
2355
2356
EXAMPLES:
2357
2358
Lets test that the del Pezzo surface `dP_6` has degree 6, as its name implies::
2359
2360
sage: dP6 = toric_varieties.dP6()
2361
sage: HH = dP6.cohomology_ring()
2362
sage: dP6.K()
2363
-V(x) - V(u) - V(y) - V(v) - V(z) - V(w)
2364
sage: dP6.integrate( HH(dP6.K())^2 )
2365
6
2366
"""
2367
from sage.schemes.toric.divisor import ToricDivisor
2368
return ToricDivisor(self, [-1]*self._fan.nrays())
2369
2370
def divisor(self, arg, base_ring=None, check=True, reduce=True):
2371
r"""
2372
Return a divisor.
2373
2374
INPUT:
2375
2376
The arguments are the same as in
2377
:func:`sage.schemes.toric.divisor.ToricDivisor`, with the
2378
exception of defining a divisor with a single integer: this method
2379
considers it to be the index of a ray of the :meth:`fan` of ``self``.
2380
2381
OUTPUT:
2382
2383
- A :class:`sage.schemes.toric.divisor.ToricDivisor_generic`
2384
2385
EXAMPLES::
2386
2387
sage: dP6 = toric_varieties.dP6()
2388
sage: dP6.coordinate_ring()
2389
Multivariate Polynomial Ring in x, u, y, v, z, w
2390
over Rational Field
2391
sage: dP6.divisor(range(6))
2392
V(u) + 2*V(y) + 3*V(v) + 4*V(z) + 5*V(w)
2393
sage: dP6.inject_variables()
2394
Defining x, u, y, v, z, w
2395
sage: dP6.divisor(x*u^3)
2396
V(x) + 3*V(u)
2397
2398
You can also construct divisors based on ray indices::
2399
2400
sage: dP6.divisor(0)
2401
V(x)
2402
sage: for i in range(0,dP6.fan().nrays()):
2403
... print dP6.divisor(i),
2404
... print ': generated by ray',
2405
... dP6.fan().ray(i)
2406
V(x) : generated by ray N(0, 1)
2407
V(u) : generated by ray N(-1, 0)
2408
V(y) : generated by ray N(-1, -1)
2409
V(v) : generated by ray N(0, -1)
2410
V(z) : generated by ray N(1, 0)
2411
V(w) : generated by ray N(1, 1)
2412
2413
TESTS:
2414
2415
We check that the issue :trac:`12812` is resolved::
2416
2417
sage: sum(dP6.divisor(i) for i in range(3))
2418
V(x) + V(u) + V(y)
2419
"""
2420
# Divisor by a ray index - must be treated here, see Trac #12812.
2421
if arg in ZZ:
2422
arg = [(1, self.gen(arg))]
2423
check = True # 1 must be coerced into the coefficient ring
2424
reduce = False
2425
from sage.schemes.toric.divisor import ToricDivisor
2426
return ToricDivisor(self, ring=base_ring, arg=arg,
2427
check=check, reduce=reduce)
2428
2429
def divisor_group(self, base_ring=ZZ):
2430
r"""
2431
Return the group of Weil divisors.
2432
2433
INPUT:
2434
2435
- ``base_ring`` -- the coefficient ring, usually ``ZZ``
2436
(default) or ``QQ``.
2437
2438
OUTPUT:
2439
2440
The (free abelian) group of Cartier divisors, that is, formal
2441
linear combinations of polynomial equations over the
2442
coefficient ring ``base_ring``.
2443
2444
These need not be toric (=defined by monomials), but allow
2445
general polynomials. The output will be an instance of
2446
:class:`sage.schemes.generic.divisor_group.DivisorGroup_generic`.
2447
2448
.. WARNING::
2449
2450
You almost certainly want the group of toric divisors, see
2451
:meth:`toric_divisor_group`. The toric divisor group is
2452
generated by the rays of the fan. The general divisor
2453
group has no toric functionality implemented.
2454
2455
EXAMPLES::
2456
2457
sage: dP6 = toric_varieties.dP6()
2458
sage: Div = dP6.divisor_group(); Div
2459
Group of ZZ-Divisors on 2-d CPR-Fano toric variety
2460
covered by 6 affine patches
2461
sage: Div(x)
2462
V(x)
2463
"""
2464
from sage.schemes.generic.divisor_group import DivisorGroup
2465
return DivisorGroup(self, base_ring)
2466
2467
def toric_divisor_group(self, base_ring=ZZ):
2468
r"""
2469
Return the group of toric (T-Weil) divisors.
2470
2471
INPUT:
2472
2473
- ``base_ring`` -- the coefficient ring, usually ``ZZ``
2474
(default) or ``QQ``.
2475
2476
OUTPUT:
2477
2478
The free Abelian agroup of toric Weil divisors, that is,
2479
formal ``base_ring``-linear combinations of codimension-one
2480
toric subvarieties. The output will be an instance of
2481
:class:`sage.schemes.toric.divisor.ToricDivisorGroup`.
2482
2483
The `i`-th generator of the divisor group is the divisor where
2484
the `i`-th homogeneous coordinate vanishes, `\{z_i=0\}`.
2485
2486
EXAMPLES::
2487
2488
sage: dP6 = toric_varieties.dP6()
2489
sage: TDiv = dP6.toric_divisor_group(); TDiv
2490
Group of toric ZZ-Weil divisors on 2-d CPR-Fano toric variety
2491
covered by 6 affine patches
2492
sage: TDiv == dP6.toric_divisor_group()
2493
True
2494
sage: TDiv.gens()
2495
(V(x), V(u), V(y), V(v), V(z), V(w))
2496
sage: dP6.coordinate_ring()
2497
Multivariate Polynomial Ring in x, u, y, v, z, w over Rational Field
2498
"""
2499
from sage.schemes.toric.divisor import ToricDivisorGroup
2500
return ToricDivisorGroup(self, base_ring);
2501
2502
def _semigroup_ring(self, cone=None, names=None):
2503
r"""
2504
Return a presentation of the semigroup ring for the dual of ``cone``.
2505
2506
INPUT:
2507
2508
See :meth:`Spec`.
2509
2510
OUTPUT:
2511
2512
For the given ``cone`` `\sigma`, return a tuple consisting of
2513
2514
* a polynomial ring `R`,
2515
2516
* an ideal `I\in R`,
2517
2518
* the dual cone `\sigma^\vee`
2519
2520
such that `R/I \sim k[\sigma^\vee \cap M]`, where `k` is the
2521
:meth:`base_ring` of the toric variety.
2522
2523
EXAMPLES::
2524
2525
sage: A2Z2 = Cone([(0,1),(2,1)])
2526
sage: AffineToricVariety(A2Z2)._semigroup_ring()
2527
(Multivariate Polynomial Ring in z0, z1, z2 over Rational Field,
2528
Ideal (-z0*z1 + z2^2) of Multivariate Polynomial Ring in z0, z1, z2 over Rational Field,
2529
2-d cone in 2-d lattice M)
2530
2531
sage: P2 = toric_varieties.P2()
2532
sage: cone = P2.fan().generating_cone(0)
2533
sage: P2._semigroup_ring(cone)
2534
(Multivariate Polynomial Ring in z0, z1 over Rational Field,
2535
Ideal (0) of Multivariate Polynomial Ring in z0, z1 over Rational Field,
2536
2-d cone in 2-d lattice M)
2537
sage: P2.change_ring(GF(101))._semigroup_ring(cone)
2538
(Multivariate Polynomial Ring in z0, z1 over Finite Field of size 101,
2539
Ideal (0) of Multivariate Polynomial Ring in z0, z1 over Finite Field of size 101,
2540
2-d cone in 2-d lattice M)
2541
"""
2542
from sage.schemes.toric.ideal import ToricIdeal
2543
if cone is None:
2544
assert self.is_affine(), \
2545
'You may only omit the cone argument for an affine toric variety!'
2546
cone = self.fan().generating_cone(0)
2547
2548
cone = self.fan().embed(cone)
2549
dual = cone.dual()
2550
basis = dual.Hilbert_basis()
2551
N = len(basis)
2552
names = normalize_names(names, N, DEFAULT_PREFIX)
2553
A = basis.column_matrix()
2554
IA = ToricIdeal(A, names, base_ring=self.base_ring())
2555
return (IA.ring(), IA, dual)
2556
2557
def Spec(self, cone=None, names=None):
2558
r"""
2559
Return the spectrum associated to the dual cone.
2560
2561
Let `\sigma \in N_\RR` be a cone and `\sigma^\vee \cap M` the
2562
associated semigroup of lattice points in the dual cone. Then
2563
2564
.. MATH::
2565
2566
S = \CC[\sigma^\vee \cap M]
2567
2568
is a `\CC`-algebra. It is spanned over `\CC` by the points of
2569
`\sigma \cap N`, addition is formal linear combination of
2570
lattice points, and multiplication of lattice points is the
2571
semigroup law (that is, addition of lattice points). The
2572
`\CC`-algebra `S` then defines a scheme `\mathop{Spec}(S)`.
2573
2574
For example, if `\sigma=\{(x,y)|x\geq 0,y\geq 0\}` is the
2575
first quadrant then `S` is the polynomial ring in two
2576
variables. The associated scheme is `\mathop{Spec}(S) =
2577
\CC^2`.
2578
2579
The same construction works over any base field, this
2580
introduction only used `\CC` for simplicity.
2581
2582
INPUT:
2583
2584
- ``cone`` -- a :class:`Cone
2585
<sage.geometry.cone.ConvexRationalPolyhedralCone>`. Can be
2586
omitted for an affine toric variety, in which case the
2587
(unique) generating cone is used.
2588
2589
- ``names`` -- (optional). Names of variables for the
2590
semigroup ring, see :func:`normalize_names` for acceptable
2591
formats. If not given, indexed variable names will be
2592
created automatically.
2593
2594
Output:
2595
2596
The spectrum of the semigroup ring `\CC[\sigma^\vee \cap M]`.
2597
2598
EXAMPLES::
2599
2600
sage: quadrant = Cone([(1,0),(0,1)])
2601
sage: AffineToricVariety(quadrant).Spec()
2602
Spectrum of Multivariate Polynomial Ring in z0, z1 over Rational Field
2603
2604
A more interesting example::
2605
2606
sage: A2Z2 = Cone([(0,1),(2,1)])
2607
sage: AffineToricVariety(A2Z2).Spec(names='u,v,t')
2608
Spectrum of Quotient of Multivariate Polynomial Ring
2609
in u, v, t over Rational Field by the ideal (-u*v + t^2)
2610
"""
2611
from sage.schemes.generic.spec import Spec
2612
R, I, dualcone = self._semigroup_ring(cone, names)
2613
return Spec(R.quotient(I))
2614
2615
def affine_algebraic_patch(self, cone=None, names=None):
2616
r"""
2617
Return the patch corresponding to ``cone`` as an affine
2618
algebraic subvariety.
2619
2620
INPUT:
2621
2622
- ``cone`` -- a :class:`Cone
2623
<sage.geometry.cone.ConvexRationalPolyhedralCone>` `\sigma`
2624
of the fan. It can be omitted for an affine toric variety,
2625
in which case the single generating cone is used.
2626
2627
OUTPUT:
2628
2629
A :class:`affine algebraic subscheme
2630
<sage.schemes.generic.algebraic_scheme.AlgebraicScheme_subscheme_affine>`
2631
corresponding to the patch `\mathop{Spec}(\sigma^\vee \cap M)`
2632
associated to the cone `\sigma`.
2633
2634
See also :meth:`affine_patch`, which expresses the patches as
2635
subvarieties of affine toric varieties instead.
2636
2637
EXAMPLES::
2638
2639
sage: cone = Cone([(0,1),(2,1)])
2640
sage: A2Z2 = AffineToricVariety(cone)
2641
sage: A2Z2.affine_algebraic_patch()
2642
Closed subscheme of Affine Space of dimension 3 over Rational Field defined by:
2643
-z0*z1 + z2^2
2644
sage: A2Z2.affine_algebraic_patch(Cone([(0,1)]), names='x, y, t')
2645
Closed subscheme of Affine Space of dimension 3 over Rational Field defined by:
2646
1
2647
"""
2648
R, I, dualcone = self._semigroup_ring(cone, names)
2649
patch_cover = AffineSpace(R)
2650
patch = patch_cover.subscheme(I)
2651
return patch
2652
2653
2654
def _orbit_closure_projection(self, cone, x):
2655
r"""
2656
Return the projection of ``x`` onto the quotient lattice of ``cone``.
2657
2658
INPUT:
2659
2660
- ``cone`` -- a :class:`cone
2661
<sage.geometry.cone.ConvexRationalPolyhedralCone>` of the :meth:`fan`
2662
of ``self``;
2663
2664
- ``x`` -- a lattice point or a cone of the :meth:`fan` of ``self``.
2665
2666
OUTPUT:
2667
2668
- the projection of ``x`` onto the quotient lattice of ``cone``, which
2669
is either a lattice point or a cone depending on the type of ``x``.
2670
This quotient lattice is the ambient lattice for the fan of the orbit
2671
closure corresponding to ``cone``.
2672
2673
If ``x`` is a cone not in the star of ``cone``, an ``IndexError`` is
2674
raised.
2675
2676
See :meth:`orbit_closure` for more details.
2677
2678
.. warning::
2679
2680
Due to incomplete support of quotient lattices (as of 12-07-2011),
2681
this function actually operates with a generic toric lattice of the
2682
same dimension as the qppropriate quotient lattice. This behaviour
2683
is likely to change in the future releases of Sage.
2684
2685
EXAMPLES::
2686
2687
sage: P2 = toric_varieties.P2()
2688
sage: H = P2.fan(1)[0]
2689
sage: [P2._orbit_closure_projection(H, p) for p in P2.fan().rays()]
2690
[(0), (1), (-1)]
2691
sage: P2._orbit_closure_projection(H, P2.fan(2)[0])
2692
1-d cone in 1-d lattice N
2693
"""
2694
cone = self.fan().embed(cone)
2695
quot = cone.sublattice_quotient()
2696
if x in cone.lattice():
2697
return vector(quot(x))
2698
2699
assert is_Cone(x)
2700
rays = [ vector(quot(r)) for r in x.rays() ]
2701
return Cone(rays)
2702
2703
# TODO: make the following work nicely.
2704
#if x in cone.lattice():
2705
#return quot(x)
2706
#assert is_Cone(x)
2707
#return Cone(x.rays(), lattice=quot)
2708
2709
2710
def orbit_closure(self, cone):
2711
r"""
2712
Return the orbit closure of ``cone``.
2713
2714
The cones `\sigma` of a fan `\Sigma` are in one-to-one correspondence
2715
with the torus orbits `O(\sigma)` of the corresponding toric variety
2716
`X_\Sigma`. Each orbit is isomorphic to a lower dimensional torus (of
2717
dimension equal to the codimension of `\sigma`). Just like the toric
2718
variety `X_\Sigma` itself, these orbits are (partially) compactified by
2719
lower-dimensional orbits. In particular, one can define the closure
2720
`V(\sigma)` of the torus orbit `O(\sigma)` in the ambient toric
2721
variety `X_\Sigma`, which is again a toric variety.
2722
2723
See Proposition 3.2.7 of [CLS]_ for more details.
2724
2725
INPUT:
2726
2727
- ``cone`` -- a :class:`cone
2728
<sage.geometry.cone.ConvexRationalPolyhedralCone>` of the fan.
2729
2730
OUTPUT:
2731
2732
- a torus orbit closure associated to ``cone`` as a
2733
:class:`toric variety <ToricVariety_field>`.
2734
2735
EXAMPLES::
2736
2737
sage: P1xP1 = toric_varieties.P1xP1()
2738
sage: H = P1xP1.fan(1)[0]
2739
sage: P1xP1.orbit_closure(H)
2740
1-d toric variety covered by 2 affine patches
2741
"""
2742
cone = self.fan().embed(cone)
2743
cones = []
2744
for star_cone in cone.star_generators():
2745
cones.append( self._orbit_closure_projection(cone, star_cone) )
2746
from sage.geometry.fan import discard_faces
2747
return ToricVariety(Fan(discard_faces(cones), check=False))
2748
2749
def count_points(self):
2750
r"""
2751
Return the number of points of ``self``.
2752
2753
This is an alias for ``point_set().cardinality()``, see
2754
:meth:`~sage.schemes.toric.homset.SchemeHomset_points_toric_field.cardinality`
2755
for details.
2756
2757
EXAMPLES::
2758
2759
sage: o = lattice_polytope.octahedron(3)
2760
sage: V = ToricVariety(FaceFan(o))
2761
sage: V2 = V.change_ring(GF(2))
2762
sage: V2.point_set().cardinality()
2763
27
2764
sage: V2.count_points()
2765
27
2766
"""
2767
return self.point_set().cardinality()
2768
2769
@cached_method
2770
def Demazure_roots(self):
2771
"""
2772
Return the Demazure roots.
2773
2774
OUTPUT:
2775
2776
The roots as points of the `M`-lattice.
2777
2778
REFERENCES:
2779
2780
.. [Demazure]
2781
M. Demazure
2782
Sous-groupes algébriques de rang maximum du groupe de Cremona.
2783
Ann. Sci. Ecole Norm. Sup. 1970, 3, 507–588.
2784
2785
.. [Bazhov]
2786
Ivan Bazhov:
2787
On orbits of the automorphism group on a complete toric variety.
2788
:arxiv:`1110.4275`,
2789
:doi:`10.1007/s13366-011-0084-0`.
2790
2791
EXAMPLE::
2792
2793
sage: P2 = toric_varieties.P2()
2794
sage: P2.Demazure_roots()
2795
(M(-1, 0), M(-1, 1), M(0, -1), M(0, 1), M(1, -1), M(1, 0))
2796
2797
Here are the remaining three examples listed in [Bazhov]_, Example 2.1 and 2.3::
2798
2799
sage: s = 3
2800
sage: cones = [(0,1),(1,2),(2,3),(3,0)]
2801
sage: Hs = ToricVariety(Fan(rays=[(1,0),(0,-1),(-1,s),(0,1)], cones=cones))
2802
sage: Hs.Demazure_roots()
2803
(M(-1, 0), M(1, 0), M(0, 1), M(1, 1), M(2, 1), M(3, 1))
2804
2805
sage: P11s = ToricVariety(Fan(rays=[(1,0),(0,-1),(-1,s)], cones=[(0,1),(1,2),(2,0)]))
2806
sage: P11s.Demazure_roots()
2807
(M(-1, 0), M(1, 0), M(0, 1), M(1, 1), M(2, 1), M(3, 1))
2808
sage: P11s.Demazure_roots() == Hs.Demazure_roots()
2809
True
2810
2811
sage: Bs = ToricVariety(Fan(rays=[(s,1),(s,-1),(-s,-1),(-s,1)], cones=cones))
2812
sage: Bs.Demazure_roots()
2813
()
2814
2815
TESTS::
2816
2817
sage: toric_varieties.A1().Demazure_roots()
2818
Traceback (most recent call last):
2819
...
2820
NotImplementedError: Demazure_roots() is only implemented for complete toric varieties.
2821
"""
2822
if not self.is_complete():
2823
raise NotImplementedError('Demazure_roots() is only implemented '
2824
'for complete toric varieties.')
2825
antiK = -self.K()
2826
fan_rays = self.fan().rays()
2827
roots = [m for m in antiK.sections()
2828
if [ray*m for ray in fan_rays].count(-1) == 1]
2829
return tuple(roots)
2830
2831
2832
def Aut_dimension(self):
2833
r"""
2834
Return the dimension of the automorphism group
2835
2836
There are three kinds of symmetries of toric varieties:
2837
2838
* Toric automorphisms (rescaling of homogeneous coordinates)
2839
2840
* Demazure roots. These are translations `x_i \to x_i +
2841
\epsilon x^m` of a homogeneous coordinate `x_i` by a
2842
monomial `x^m` of the same homogeneous degree.
2843
2844
* Symmetries of the fan. These yield discrete subgroups.
2845
2846
OUTPUT:
2847
2848
An integer. The dimension of the automorphism group. Equals
2849
the dimension of the `M`-lattice plus the number of Demazure
2850
roots.
2851
2852
EXAMPLES::
2853
2854
sage: P2 = toric_varieties.P2()
2855
sage: P2.Aut_dimension()
2856
8
2857
2858
TESTS::
2859
2860
sage: toric_varieties.A1().Aut_dimension()
2861
Traceback (most recent call last):
2862
...
2863
NotImplementedError: Aut_dimension() is only implemented for complete toric varieties.
2864
"""
2865
if not self.is_complete():
2866
raise NotImplementedError('Aut_dimension() is only implemented '
2867
'for complete toric varieties.')
2868
return self.fan().lattice_dim() + len(self.Demazure_roots())
2869
2870
2871
def normalize_names(names=None, ngens=None, prefix=None, indices=None,
2872
return_prefix=False):
2873
r"""
2874
Return a list of names in the standard form.
2875
2876
INPUT:
2877
2878
All input parameters are optional.
2879
2880
- ``names`` -- names given either as a single string (with individual
2881
names separated by commas or spaces) or a list of strings with each
2882
string specifying a name. If the last name ends with the plus sign,
2883
"+", this name will be used as ``prefix`` (even if ``prefix`` was
2884
given explicitly);
2885
2886
- ``ngens`` -- number of names to be returned;
2887
2888
- ``prefix`` -- prefix for the indexed names given as a string;
2889
2890
- ``indices`` -- list of integers (default: ``range(ngens)``) used as
2891
indices for names with ``prefix``. If given, must be of length
2892
``ngens``;
2893
2894
- ``return_prefix`` -- if ``True``, the last element of the returned list
2895
will contain the prefix determined from ``names`` or given as the
2896
parameter ``prefix``. This is useful if you may need more names in the
2897
future.
2898
2899
OUTPUT:
2900
2901
- list of names given as strings.
2902
2903
These names are constructed in the following way:
2904
2905
#. If necessary, split ``names`` into separate names.
2906
#. If the last name ends with "+", put it into ``prefix``.
2907
#. If ``ngens`` was given, add to the names obtained so far as many
2908
indexed names as necessary to get this number. If the ``k``-th name of
2909
the *total* list of names is indexed, it is
2910
``prefix + str(indices[k])``. If there were already more names than
2911
``ngens``, discard "extra" ones.
2912
#. Check if constructed names are valid. See :func:`certify_names` for
2913
details.
2914
#. If the option ``return_prefix=True`` was given, add ``prefix`` to the
2915
end of the list.
2916
2917
EXAMPLES:
2918
2919
As promised, all parameters are optional::
2920
2921
sage: from sage.schemes.toric.variety import normalize_names
2922
sage: normalize_names()
2923
[]
2924
2925
One of the most common uses is probably this one::
2926
2927
sage: normalize_names("x+", 4)
2928
['x0', 'x1', 'x2', 'x3']
2929
2930
Now suppose that you want to enumerate your variables starting with one
2931
instead of zero::
2932
2933
sage: normalize_names("x+", 4, indices=range(1,5))
2934
['x1', 'x2', 'x3', 'x4']
2935
2936
You may actually have an arbitrary enumeration scheme::
2937
2938
sage: normalize_names("x+", 4, indices=[1, 10, 100, 1000])
2939
['x1', 'x10', 'x100', 'x1000']
2940
2941
Now let's add some "explicit" names::
2942
2943
sage: normalize_names("x y z t+", 4)
2944
['x', 'y', 'z', 't3']
2945
2946
Note that the "automatic" name is ``t3`` instead of ``t0``. This may seem
2947
weird, but the reason for this behaviour is that the fourth name in this
2948
list will be the same no matter how many explicit names were given::
2949
2950
sage: normalize_names("x y t+", 4)
2951
['x', 'y', 't2', 't3']
2952
2953
This is especially useful if you get ``names`` from a user but want to
2954
specify all default names::
2955
2956
sage: normalize_names("x, y", 4, prefix="t")
2957
['x', 'y', 't2', 't3']
2958
2959
In this format, the user can easily override your choice for automatic
2960
names::
2961
2962
sage: normalize_names("x y s+", 4, prefix="t")
2963
['x', 'y', 's2', 's3']
2964
2965
Let's now use all parameters at once::
2966
2967
sage: normalize_names("x, y, s+", 4, prefix="t",
2968
... indices=range(1,5), return_prefix=True)
2969
['x', 'y', 's3', 's4', 's']
2970
2971
Note that you still need to give indices for all names, even if some of
2972
the first ones will be "wasted" because of the explicit names. The reason
2973
is the same as before - this ensures consistency of automatically
2974
generated names, no matter how many explicit names were given.
2975
2976
The prefix is discarded if ``ngens`` was not given::
2977
2978
sage: normalize_names("alpha, beta, gamma, zeta+")
2979
['alpha', 'beta', 'gamma']
2980
2981
Finally, let's take a look at some possible mistakes::
2982
2983
sage: normalize_names("123")
2984
Traceback (most recent call last):
2985
...
2986
ValueError: name must start with a letter! Got 123
2987
2988
A more subtle one::
2989
2990
sage: normalize_names("x1", 4, prefix="x")
2991
Traceback (most recent call last):
2992
...
2993
ValueError: names must be distinct! Got: ['x1', 'x1', 'x2', 'x3']
2994
"""
2995
if names is None:
2996
names = []
2997
elif isinstance(names, str):
2998
names = names.replace(",", " ").split()
2999
else:
3000
try:
3001
names = list(names)
3002
except TypeError:
3003
raise TypeError(
3004
"names must be a string or a list or tuple of them!")
3005
for name in names:
3006
if not isinstance(name, str):
3007
raise TypeError(
3008
"names must be a string or a list or tuple of them!")
3009
if names and names[-1].endswith("+"):
3010
prefix = names.pop()[:-1]
3011
if ngens is None:
3012
ngens = len(names)
3013
if len(names) < ngens:
3014
if prefix is None:
3015
raise IndexError("need %d names but only %d are given!"
3016
% (ngens, len(names)))
3017
if indices is None:
3018
indices = range(ngens)
3019
elif len(indices) != ngens:
3020
raise ValueError("need exactly %d indices, but got %d!"
3021
% (ngens, len(indices)))
3022
names += [prefix + str(i) for i in indices[len(names):]]
3023
if len(names) > ngens:
3024
names = names[:ngens]
3025
# Check that all given and constructed names are valid
3026
certify_names(names)
3027
if return_prefix:
3028
names.append(prefix)
3029
return names
3030
3031
3032
def certify_names(names):
3033
r"""
3034
Make sure that ``names`` are valid in Python.
3035
3036
INPUT:
3037
3038
- ``names`` -- list of strings.
3039
3040
OUTPUT:
3041
3042
- none, but a ``ValueError`` exception is raised if ``names`` are invalid.
3043
3044
Each name must satisfy the following requirements:
3045
3046
* Be non-empty.
3047
* Contain only (Latin) letters, digits, and underscores ("_").
3048
* Start with a letter.
3049
3050
In addition, all names must be distinct.
3051
3052
EXAMPLES::
3053
3054
sage: from sage.schemes.toric.variety import certify_names
3055
sage: certify_names([])
3056
sage: certify_names(["a", "x0", "x_45"])
3057
sage: certify_names(["", "x0", "x_45"])
3058
Traceback (most recent call last):
3059
...
3060
ValueError: name must be nonempty!
3061
sage: certify_names(["a", "0", "x_45"])
3062
Traceback (most recent call last):
3063
...
3064
ValueError: name must start with a letter! Got 0
3065
sage: certify_names(["a", "x0", "@_45"])
3066
Traceback (most recent call last):
3067
...
3068
ValueError: name must be alphanumeric! Got @_45
3069
sage: certify_names(["a", "x0", "x0"])
3070
Traceback (most recent call last):
3071
...
3072
ValueError: names must be distinct! Got: ['a', 'x0', 'x0']
3073
"""
3074
for name in names:
3075
if not name:
3076
raise ValueError("name must be nonempty!")
3077
if not name.isalnum() and not name.replace("_","").isalnum():
3078
# Must be alphanumeric except for non-leading '_'
3079
raise ValueError("name must be alphanumeric! Got %s" % name)
3080
if not name[0].isalpha():
3081
raise ValueError("name must start with a letter! Got %s" % name)
3082
if len(set(names)) != len(names):
3083
raise ValueError("names must be distinct! Got: %s" % names)
3084
3085
3086
#*****************************************************************
3087
class CohomologyRing(QuotientRing_generic, UniqueRepresentation):
3088
r"""
3089
The (even) cohomology ring of a toric variety.
3090
3091
Irregardles of the variety's base ring, we always work with the
3092
variety over `\CC` and its topology.
3093
3094
The cohomology is always the singular cohomology with
3095
`\QQ`-coefficients. Note, however, that the cohomology of smooth
3096
toric varieties is torsion-free, so there is no loss of
3097
information in that case.
3098
3099
Currently, the toric variety must not be "too singular". See
3100
:meth:`ToricVariety_field.cohomology_ring` for a detailed
3101
description of which toric varieties are admissible. For such
3102
varieties the odd-dimensional cohomology groups vanish.
3103
3104
.. WARNING::
3105
3106
You should not create instances of this class manually. Use
3107
:meth:`ToricVariety_field.cohomology_ring` to generate the
3108
cohomology ring.
3109
3110
INPUT:
3111
3112
- ``variety`` -- a toric variety. Currently, the toric variety
3113
must be at least an orbifold. See
3114
:meth:`ToricVariety_field.cohomology_ring` for a detailed
3115
description of which toric varieties are admissible.
3116
3117
EXAMPLES::
3118
3119
sage: P2 = toric_varieties.P2()
3120
sage: P2.cohomology_ring()
3121
Rational cohomology ring of a 2-d CPR-Fano toric variety covered by 3 affine patches
3122
3123
This is equivalent to::
3124
3125
sage: from sage.schemes.toric.variety import CohomologyRing
3126
sage: CohomologyRing(P2)
3127
Rational cohomology ring of a 2-d CPR-Fano toric variety covered by 3 affine patches
3128
"""
3129
3130
def __init__(self, variety):
3131
r"""
3132
See :class:`CohomologyRing` for documentation.
3133
3134
TESTS::
3135
3136
sage: P2 = toric_varieties.P2()
3137
sage: P2.cohomology_ring()
3138
Rational cohomology ring of a 2-d CPR-Fano toric variety covered by 3 affine patches
3139
3140
TESTS::
3141
3142
sage: cone1 = Cone([(1,0)]); cone2 = Cone([(1,0)])
3143
sage: cone1 is cone2
3144
False
3145
sage: fan1 = Fan([cone1]); fan2 = Fan([cone2])
3146
sage: fan1 is fan2
3147
False
3148
sage: X1 = ToricVariety(fan1); X2 = ToricVariety(fan2)
3149
sage: X1 is X2
3150
False
3151
sage: X1.cohomology_ring() is X2.cohomology_ring() # see http://trac.sagemath.org/sage_trac/ticket/10325
3152
True
3153
sage: TDiv = X1.toric_divisor_group()
3154
sage: X1.toric_divisor_group() is TDiv
3155
True
3156
sage: X2.toric_divisor_group() is TDiv
3157
True
3158
sage: TDiv.scheme() is X1 # as you expect
3159
True
3160
sage: TDiv.scheme() is X2 # perhaps less obvious, but toric_divisor_group is unique!
3161
False
3162
sage: TDiv.scheme() == X2 # isomorphic, but not necessarily identical
3163
True
3164
sage: TDiv.scheme().cohomology_ring() is X2.cohomology_ring() # this is where it gets tricky
3165
True
3166
sage: TDiv.gen(0).Chern_character() * X2.cohomology_ring().one()
3167
[1]
3168
"""
3169
self._variety = variety
3170
3171
if not variety.is_orbifold():
3172
raise NotImplementedError, 'Requires an orbifold toric variety.'
3173
3174
R = PolynomialRing(QQ, variety.variable_names())
3175
self._polynomial_ring = R
3176
3177
I = variety._fan.linear_equivalence_ideal(R) + variety._fan.Stanley_Reisner_ideal(R)
3178
super(CohomologyRing, self).__init__(R, I, names=variety.variable_names())
3179
3180
def _repr_(self):
3181
r"""
3182
Return a string representation of the cohomology ring.
3183
3184
OUTPUT:
3185
3186
A string.
3187
3188
EXAMPLES::
3189
3190
sage: toric_varieties.P2().cohomology_ring()._repr_()
3191
'Rational cohomology ring of a 2-d CPR-Fano toric variety covered by 3 affine patches'
3192
"""
3193
return 'Rational cohomology ring of a '+self._variety._repr_()
3194
3195
def _latex_(self):
3196
r"""
3197
Return a latex representation of the cohomology ring.
3198
3199
OUTPUT:
3200
3201
A string.
3202
3203
EXAMPLES::
3204
3205
sage: cohomology_ring = toric_varieties.P2().cohomology_ring()
3206
sage: cohomology_ring._latex_()
3207
'H^\\ast\\left(\\mathbb{P}_{\\Delta^{2}},\\QQ\\right)'
3208
"""
3209
return 'H^\\ast\\left('+self._variety._latex_()+',\QQ\\right)'
3210
3211
def _element_constructor_(self,x):
3212
r"""
3213
Construct a :class:`CohomologyClass`.
3214
3215
INPUT::
3216
3217
- ``x`` -- something that defines a cohomology class. Either a
3218
cohomology class, a cone of the fan, or something that can
3219
be converted into a polynomial in the homogeneous
3220
coordinates.
3221
3222
OUTPUT:
3223
3224
The :class:`CohomologyClass` defined by ``x``.
3225
3226
EXAMPLES::
3227
3228
sage: dP6 = toric_varieties.dP6()
3229
sage: H = dP6.cohomology_ring()
3230
sage: cone = dP6.fan().cone_containing(2,3); cone
3231
2-d cone of Rational polyhedral fan in 2-d lattice N
3232
sage: H(cone) # indirect doctest
3233
[-w^2]
3234
sage: H( Cone(cone) )
3235
[-w^2]
3236
sage: H( dP6.fan(0)[0] ) # trivial cone
3237
[1]
3238
3239
Non-smooth cones are a bit tricky. For such a cone, the
3240
intersection of the divisors corresponding to the rays is
3241
still proportional to the product of the variables, but the
3242
coefficient is a multiple depending on the orbifold
3243
singularity. See also [CLS]_, Lemma 12.5.2::
3244
3245
sage: P2_123 = toric_varieties.P2_123()
3246
sage: HH = P2_123.cohomology_ring()
3247
sage: HH(Cone([(1,0)])) * HH(Cone([(-2,-3)]))
3248
[2*z2^2]
3249
sage: HH(Cone([(1,0), (-2,-3)]))
3250
[6*z2^2]
3251
sage: [ HH(c) for c in P2_123.fan().generating_cones() ]
3252
[[6*z2^2], [6*z2^2], [6*z2^2]]
3253
sage: P2_123.volume_class()
3254
[6*z2^2]
3255
sage: [ HH(c.facets()[0]) * HH(c.facets()[1]) for c in P2_123.fan().generating_cones() ]
3256
[[6*z2^2], [3*z2^2], [2*z2^2]]
3257
3258
Numbers will be converted into the ring::
3259
3260
sage: P2 = toric_varieties.P2()
3261
sage: H = P2.cohomology_ring()
3262
sage: H._element_constructor_(1)
3263
[1]
3264
sage: H(1)
3265
[1]
3266
sage: type( H(1) )
3267
<class 'sage.schemes.toric.variety.CohomologyClass'>
3268
sage: P2.inject_variables()
3269
Defining x, y, z
3270
sage: H(1+x*y+z)
3271
[z^2 + z + 1]
3272
"""
3273
fan = self._variety.fan()
3274
if isinstance(x, CohomologyClass) and x.parent()==self:
3275
return x
3276
if isinstance(x, QuotientRingElement):
3277
x = x.lift()
3278
elif is_Cone(x):
3279
cone = fan.embed(x)
3280
assert cone.ambient() is fan
3281
mult = cone.rays().column_matrix().index_in_saturation()
3282
x = prod((self.cover_ring().gen(i) for i in cone.ambient_ray_indices()),
3283
z=self.cover_ring().one()) * mult
3284
else:
3285
try:
3286
# divisor, for example, know how to compute their own cohomology class
3287
return x.cohomology_class()
3288
except AttributeError:
3289
# this ensures that rationals are converted to cohomology ring elements
3290
x = self.cover_ring()(x)
3291
return CohomologyClass(self, x)
3292
3293
# We definitely should not override __call__, but since our
3294
# superclass QuotientRing_generic does not adhere to the coercion
3295
# model we cannot either. See
3296
# http://trac.sagemath.org/sage_trac/ticket/9429
3297
def __call__(self, x, coerce=True):
3298
r"""
3299
Turn ``x`` into a ``CohomologyClass``.
3300
3301
EXAMPLES::
3302
3303
sage: P2 = toric_varieties.P2()
3304
sage: H = P2.cohomology_ring()
3305
sage: H(1)
3306
[1]
3307
sage: type( H(1) )
3308
<class 'sage.schemes.toric.variety.CohomologyClass'>
3309
"""
3310
return self._element_constructor_(x)
3311
3312
def gens(self):
3313
r"""
3314
Return the generators of the cohomology ring.
3315
3316
OUTPUT:
3317
3318
A tuple of generators, one for each toric divisor of the toric
3319
variety ``X``. The order is the same as the ordering of the
3320
rays of the fan ``X.fan().rays()``, which is also the same as
3321
the ordering of the one-cones in ``X.fan(1)``
3322
3323
EXAMPLES::
3324
3325
sage: P2 = toric_varieties.P2()
3326
sage: P2.cohomology_ring().gens()
3327
([z], [z], [z])
3328
"""
3329
if "_gens" not in self.__dict__:
3330
self._gens = tuple( self.gen(i) for i in range(0,self._variety.fan().nrays()) )
3331
return self._gens
3332
3333
def gen(self, i):
3334
r"""
3335
Return the generators of the cohomology ring.
3336
3337
INPUT:
3338
3339
- ``i`` -- integer.
3340
3341
OUTPUT:
3342
3343
The ``i``-th generator of the cohomology ring. If we denote
3344
the toric variety by ``X``, then this generator is
3345
associated to the ray ``X.fan().ray(i)``, which spans the
3346
one-cone ``X.fan(1)[i]``
3347
3348
EXAMPLES::
3349
3350
sage: P2 = toric_varieties.P2()
3351
sage: P2.cohomology_ring().gen(2)
3352
[z]
3353
"""
3354
return CohomologyClass(self, self._polynomial_ring.gen(i))
3355
3356
3357
#*****************************************************************
3358
def is_CohomologyClass(x):
3359
r"""
3360
Check whether ``x`` is a cohomology class of a toric variety.
3361
3362
INPUT:
3363
3364
- ``x`` -- anything.
3365
3366
OUTPUT:
3367
3368
``True`` or ``False`` depending on whether ``x`` is an instance of
3369
:class:`CohomologyClass`
3370
3371
EXAMPLES::
3372
3373
sage: P2 = toric_varieties.P2()
3374
sage: HH = P2.cohomology_ring()
3375
sage: from sage.schemes.toric.variety import is_CohomologyClass
3376
sage: is_CohomologyClass( HH.one() )
3377
True
3378
sage: is_CohomologyClass( HH(P2.fan(1)[0]) )
3379
True
3380
sage: is_CohomologyClass('z')
3381
False
3382
"""
3383
return isinstance(x,CohomologyClass)
3384
3385
3386
#*****************************************************************
3387
class CohomologyClass(QuotientRingElement):
3388
r"""
3389
An element of the :class:`CohomologyRing`.
3390
3391
.. WARNING::
3392
3393
You should not create instances of this class manually. The
3394
generators of the cohomology ring as well as the cohomology
3395
classes associated to cones of the fan can be obtained from
3396
:meth:`ToricVariety_field.cohomology_ring`.
3397
3398
EXAMPLES::
3399
3400
sage: P2 = toric_varieties.P2()
3401
sage: P2.cohomology_ring().gen(0)
3402
[z]
3403
sage: HH = P2.cohomology_ring()
3404
sage: HH.gen(0)
3405
[z]
3406
sage: cone = P2.fan(1)[0]; HH(cone)
3407
[z]
3408
"""
3409
3410
def __init__(self, cohomology_ring, representative):
3411
r"""
3412
Construct the cohomology class.
3413
3414
INPUT:
3415
3416
- ``cohomology_ring`` -- :class:`CohomologyRing`.
3417
3418
- ``representative`` -- a polynomial in the generators of the cohomology ring.
3419
3420
OUTPUT:
3421
3422
An instance of :class:`CohomologyClass`.
3423
3424
EXAMPLES::
3425
3426
sage: P2 = toric_varieties.P2()
3427
sage: H = P2.cohomology_ring()
3428
sage: from sage.schemes.toric.variety import CohomologyClass
3429
sage: CohomologyClass(H, H.defining_ideal().ring().zero() )
3430
[0]
3431
"""
3432
assert representative in cohomology_ring.defining_ideal().ring(), \
3433
'The given representative is not in the parent polynomial ring.'
3434
super(CohomologyClass, self).__init__(cohomology_ring, representative)
3435
3436
def _repr_(self):
3437
r"""
3438
Return a string representation of the cohomology class.
3439
3440
OUTPUT:
3441
3442
A string.
3443
3444
EXAMPLES::
3445
3446
sage: toric_varieties.P2().cohomology_ring().gen(0)._repr_()
3447
'[z]'
3448
"""
3449
return '['+super(CohomologyClass,self)._repr_()+']'
3450
3451
def _latex_(self):
3452
r"""
3453
Return a latex representation of the cohomology class.
3454
3455
OUTPUT:
3456
3457
A string.
3458
3459
EXAMPLES::
3460
3461
sage: cohomology_class = toric_varieties.P2().cohomology_ring().gen(0)^2/2
3462
sage: cohomology_class._latex_()
3463
'\\left[ \\frac{1}{2} z^{2} \\right]'
3464
"""
3465
return r'\left[ %s \right]' % latex(self.lift())
3466
3467
def deg(self):
3468
r"""
3469
The degree of the cohomology class.
3470
3471
OUTPUT:
3472
3473
An integer `d` such that the cohomology class is in degree
3474
`2d`. If the cohomology class is of mixed degree, the highest
3475
degree is returned.
3476
3477
EXAMPLES::
3478
3479
sage: P2 = toric_varieties.P2()
3480
sage: P2.cohomology_ring().gen(0).deg()
3481
1
3482
sage: P2.cohomology_ring().zero().deg()
3483
-1
3484
"""
3485
return self.lift().degree()
3486
3487
def part_of_degree(self, d):
3488
r"""
3489
Project the (mixed-degree) cohomology class to the given degree.
3490
3491
.. MATH::
3492
3493
\mathop{pr}\nolimits_d:~ H^\bullet(X_\Delta,\QQ) \to H^{2d}(X_\Delta,\QQ)
3494
3495
INPUT:
3496
3497
- An integer ``d``
3498
3499
OUTPUT:
3500
3501
- The degree-``2d`` part of the cohomology class.
3502
3503
EXAMPLES::
3504
3505
sage: P1xP1 = toric_varieties.P1xP1()
3506
sage: t = P1xP1.cohomology_ring().gen(0)
3507
sage: y = P1xP1.cohomology_ring().gen(2)
3508
sage: 3*t+4*t^2*y+y+t*y+t+1
3509
[t*y + 4*t + y + 1]
3510
sage: (3*t+4*t^2*y+y+t*y+t+1).part_of_degree(1)
3511
[4*t + y]
3512
"""
3513
Q = self.parent()
3514
# We iterate over monomials of self.lift()
3515
p = filter( lambda x: x[1].total_degree() == d, self.lift() )
3516
if len(p)==0:
3517
return Q.zero()
3518
else:
3519
return Q(sum(x[0]*x[1] for x in p))
3520
3521
def exp(self):
3522
"""
3523
Exponentiate ``self``.
3524
3525
.. NOTE::
3526
3527
The exponential `\exp(x)` of a rational number `x` is
3528
usually not rational. Therefore, the cohomology class must
3529
not have a constant (degree zero) part. The coefficients
3530
in the Taylor series of `\exp` are rational, so any
3531
cohomology class without constant term can be
3532
exponentiated.
3533
3534
OUTPUT
3535
3536
The cohomology class `\exp(` ``self`` `)` if the constant part
3537
vanishes, otherwise a ``ValueError`` is raised.
3538
3539
EXAMPLES::
3540
3541
sage: P2 = toric_varieties.P2()
3542
sage: H_class = P2.cohomology_ring().gen(0)
3543
sage: H_class
3544
[z]
3545
sage: H_class.exp()
3546
[1/2*z^2 + z + 1]
3547
"""
3548
if not self.part_of_degree(0).is_zero():
3549
raise ValueError, 'Must not have a constant part.'
3550
exp_x = self.parent().one()
3551
for d in range(1,self.parent()._variety.dimension()+1):
3552
exp_x += self**d / factorial(d)
3553
return exp_x
3554
3555
3556