Macaulay2 » Documentation
Packages » NumericalImplicitization :: pseudoWitnessSet
next | previous | forward | backward | up | index | toc

pseudoWitnessSet -- computes a pseudo-witness set for the image of a variety

Synopsis

Description

This method computes a pseudo-witness set for the image of a variety, by computing the intersection of the image with a complementary-dimensional linear slice via tracking monodromy loops with homotopy continuation, and then applying the trace test. If the trace test fails, only a lower bound for the degree and an incomplete pseudo-witness set is returned. This technique circumvents the calculation of the kernel of the associated ring map.

The method also allows the user to provide a particular linear slice $L$ of the image. In this case a list of point pairs $(p, q)$ such that $p$ is in $V(I)$, $q = F(p)$, and $q$ is in $L$, must be provided (to have an initial input point to the monodromy - even if it only consists of a single such pair). The method then applies monodromy to try to compute the entire intersection $F(V(I))\cap L$. If no linear slice is given, then a random complementary-dimensional linear slice will be chosen, in which case no seed is needed, as an initial point pair will be chosen to lie on the slice.

The following example computes the degree of the Grassmannian $Gr(2,4)$ of $P^1$'s in $P^3$, under its Plücker embedding in $P^5$.

i1 : R = CC[x_(1,1)..x_(2,4)];
i2 : F = (minors(2, genericMatrix(R, 2, 4)))_*;
i3 : W = pseudoWitnessSet(F, ideal 0_R)
Sampling point in source ...
Tracking monodromy loops ...
Points found: 2
Points found: 2
Points found: 2
Points found: 2
Running trace test ...

o3 = a "pseudo-witness set", indicating
     the image has degree 2

o3 : PseudoWitnessSet
i4 : W.isCompletePseudoWitnessSet

o4 = true
i5 : W.degree

o5 = 2

This method can also handle cases where the parameterization has positive dimensional fibers. In the example below, we verify that the variety of $3 x 3 x 3$ tensors of border rank $<= 4$, i.e. the $4$th secant variety of $P^2 x P^2 x P^2$, has degree $9$. This is a hypersurface, with defining equation known as Strassen's invariant, and it is also a defective secant variety (meaning its dimension is less than expected). Here, the parametrization has $10$ dimensional fibers. For more on this example, see V. Strassen, $The asymptotic spectrum of tensors$, J. Reine Angew. Math. 384 (1988), 102-152.

i6 : R = CC[a_(0,0)..a_(3,2), b_(0,0)..b_(3,2), c_(0,0)..c_(3,2)];
i7 : F = toList apply((0,0,0)..(2,2,2), (i,j,k) ->
        a_(0,i)*b_(0,j)*c_(0,k) +
        a_(1,i)*b_(1,j)*c_(1,k) +
        a_(2,i)*b_(2,j)*c_(2,k) +
        a_(3,i)*b_(3,j)*c_(3,k));
i8 : pseudoWitnessSet(F, ideal 0_R, Repeats => 2)
Sampling point in source ...
Tracking monodromy loops ...
Points found: 1
Points found: 2
Points found: 3
Points found: 5
Points found: 7
Points found: 9
Points found: 9
Points found: 9
Running trace test ...

o8 = a pseudo-witness set, indicating
    the degree of the image is 9

o8 : PseudoWitnessSet

Finally, this method has a large number of optional inputs which may be specified by the user to fit a particular problem instance.

The option Repeats sets the maximum number of consecutive repetitive monodromy loops when computing a pseudo-witness set. A repetitive monodromy loop is one where no new points in the image are discovered. After this many consecutive repetitive monodromy loops occur, the trace test is applied to determine if a complete pseudo-witness set has been found. The default value is $3$.

The option MaxAttempts sets the maximum number of times the trace test will be attempted when computing a pseudo-witness set. After a trace test fails, a new slice is chosen, the previous points are tracked to the new slice, and monodromy is performed anew. If the trace test has failed MaxAttempts many times, an incomplete pseudo-witness set is returned. The default value is $5$.

Here is an example in which a badly chosen random seed results in a failed trace test on the first attempt. In later attempts, the trace test passes and the degree of the twisted cubic is correctly computed to be $3$.

i6 : setRandomSeed 10

o6 = 10
i7 : R = CC[s,t]

o7 = R

o7 : PolynomialRing
i8 : F = basis(3, R)
-- warning: experimental computation over inexact field begun
--          results not reliable (one warning given per session)

o8 = | s3 s2t st2 t3 |

             1      4
o8 : Matrix R  <-- R
i9 : pseudoWitnessSet(F, ideal 0_R)
Sampling point in source ...
Tracking monodromy loops ...
Points found: 2
Points found: 2
Points found: 2
Points found: 2
Running trace test ...
Failed trace test! Trace: .0128189
Tracking monodromy loops ...
Points found: 2
Points found: 3
Points found: 3
Points found: 3
Points found: 3
Running trace test ...

o9 = a "pseudo-witness set", indicating
     the image has degree 3

o9 : PseudoWitnessSet

We compare this with the native $Macaulay2$ function degree (using a symbolic Gr&ouml;bner basis computation).

i10 : degree ker map(QQ[s,t], QQ[y_0..y_3], {s^3,s^2*t,s*t^2,t^3})

o10 = 3

The option MaxPoints sets a number of points such that if more than this number of points is found following a monodromy loop, then the method gracefully exits. The option is especially useful in the case that the user specifies a linear slice $L$ (as discussed above) which is in special position with respect to $F(V(I))$ (e.g. if $F(V(I))\cap L$ is positive-dimensional). The default value is infinity.

The option DoRefinements specifies whether or not to refine solution points found via monodromy. Refinement of points may improve their accuracy. If the value of this option is true, then refinement occurs after every tracking (which may increase the time for computation). The default value is false.

The option DoTraceTest specifies whether or not to run the trace test. This is useful when the user specifies a special linear slice $L$ (as in the discussion on MaxPoints above). The default value is true.

The option TraceThreshold sets the threshold for a pseudo-witness set to pass the trace test. The trace for a complete exact pseudo-witness set is $0$; large nonzero values indicate failure (the larger the value, the worse the failure). The default value is $1e-5$.

The option Threshold sets the threshold for determining point equality. If this option has value $n$, then two points are considered equal iff their first $n$ significant digits agree (equivalently, in scientific notation, the exponents and first $n$ digits of the mantissa agree). The default value is $5$.

See also

Ways to use pseudoWitnessSet :

For the programmer

The object pseudoWitnessSet is a method function with options.