001/*
002 * $Id$
003 */
004
005package edu.jas.root;
006
007
008import java.util.ArrayList;
009import java.util.List;
010import java.util.Map;
011import java.util.Iterator;
012import java.util.TreeMap;
013import java.util.SortedMap;
014
015import org.apache.logging.log4j.Logger;
016import org.apache.logging.log4j.LogManager; 
017
018import edu.jas.arith.BigDecimal;
019import edu.jas.arith.BigRational;
020import edu.jas.arith.Rational;
021import edu.jas.poly.Complex;
022import edu.jas.poly.ComplexRing;
023import edu.jas.poly.GenPolynomial;
024import edu.jas.poly.GenPolynomialRing;
025import edu.jas.poly.PolyUtil;
026import edu.jas.structure.RingElem;
027import edu.jas.structure.RingFactory;
028import edu.jas.structure.UnaryFunctor;
029import edu.jas.ufd.Squarefree;
030import edu.jas.ufd.SquarefreeFactory;
031
032
033/**
034 * Complex roots abstract class.
035 * @param <C> coefficient type.
036 * @author Heinz Kredel
037 */
038public abstract class ComplexRootsAbstract<C extends RingElem<C> & Rational> implements ComplexRoots<C> {
039
040
041    private static final Logger logger = LogManager.getLogger(ComplexRootsAbstract.class);
042
043
044    private static final boolean debug = logger.isDebugEnabled();
045
046
047    /**
048     * Engine for square free decomposition.
049     */
050    public final Squarefree<Complex<C>> engine;
051
052
053    /**
054     * Constructor.
055     * @param cf coefficient factory.
056     */
057    public ComplexRootsAbstract(RingFactory<Complex<C>> cf) {
058        if (!(cf instanceof ComplexRing)) {
059            throw new IllegalArgumentException("cf not supported coefficients " + cf);
060        }
061        engine = SquarefreeFactory.<Complex<C>> getImplementation(cf);
062    }
063
064
065    /**
066     * Root bound. With f(-M + i M) * f(-M - i M) * f(M - i M) * f(M + i M) !=
067     * 0.
068     * @param f univariate polynomial.
069     * @return M such that root(f) is contained in the rectangle spanned by M.
070     */
071    public Complex<C> rootBound(GenPolynomial<Complex<C>> f) {
072        if (f == null) {
073            return null;
074        }
075        RingFactory<Complex<C>> cfac = f.ring.coFac;
076        Complex<C> M = cfac.getONE();
077        if (f.isZERO() || f.isConstant()) {
078            return M;
079        }
080        Complex<C> a = f.leadingBaseCoefficient().norm();
081        for (Complex<C> c : f.getMap().values()) {
082            Complex<C> d = c.norm().divide(a);
083            if (M.compareTo(d) < 0) {
084                M = d;
085            }
086        }
087        M = M.sum(cfac.getONE());
088        //System.out.println("M = " + M);
089        return M;
090    }
091
092
093    /**
094     * Magnitude bound.
095     * @param rect rectangle.
096     * @param f univariate polynomial.
097     * @return B such that |f(c)| &lt; B for c in rect.
098     */
099    public C magnitudeBound(Rectangle<C> rect, GenPolynomial<Complex<C>> f) {
100        if (f == null) {
101            return null;
102        }
103        if (f.isZERO()) {
104            return f.ring.coFac.getONE().getRe();
105        }
106        //System.out.println("f = " + f);
107        if (f.isConstant()) {
108            Complex<C> c = f.leadingBaseCoefficient();
109            return c.norm().getRe();
110        }
111        GenPolynomial<Complex<C>> fa = f.map(new UnaryFunctor<Complex<C>, Complex<C>>() {
112
113
114                public Complex<C> eval(Complex<C> a) {
115                    return a.norm();
116                }
117            });
118        //System.out.println("fa = " + fa);
119        Complex<C> Mc = rect.getNW().norm();
120        C M = Mc.getRe();
121        //System.out.println("M = " + M);
122        Complex<C> M1c = rect.getSW().norm();
123        C M1 = M1c.getRe();
124        if (M.compareTo(M1) < 0) {
125            M = M1;
126            Mc = M1c;
127        }
128        M1c = rect.getSE().norm();
129        M1 = M1c.getRe();
130        if (M.compareTo(M1) < 0) {
131            M = M1;
132            Mc = M1c;
133        }
134        M1c = rect.getNE().norm();
135        M1 = M1c.getRe();
136        if (M.compareTo(M1) < 0) {
137            //M = M1;
138            Mc = M1c;
139        }
140        //System.out.println("M = " + M);
141        Complex<C> B = PolyUtil.<Complex<C>> evaluateMain(f.ring.coFac, fa, Mc);
142        //System.out.println("B = " + B);
143        return B.getRe();
144    }
145
146
147    /**
148     * Complex root count of complex polynomial on rectangle.
149     * @param rect rectangle.
150     * @param a univariate complex polynomial.
151     * @return root count of a in rectangle.
152     */
153    public abstract long complexRootCount(Rectangle<C> rect, GenPolynomial<Complex<C>> a)
154        throws InvalidBoundaryException;
155
156
157    /**
158     * List of complex roots of complex polynomial a on rectangle.
159     * @param rect rectangle.
160     * @param a univariate squarefree complex polynomial.
161     * @return list of complex roots.
162     */
163    public abstract List<Rectangle<C>> complexRoots(Rectangle<C> rect, GenPolynomial<Complex<C>> a)
164        throws InvalidBoundaryException;
165
166
167    /**
168     * List of complex roots of complex polynomial.
169     * @param a univariate complex polynomial.
170     * @return list of complex roots.
171     */
172    @SuppressWarnings("unchecked")
173    public List<Rectangle<C>> complexRoots(GenPolynomial<Complex<C>> a) {
174        List<Rectangle<C>> roots = new ArrayList<Rectangle<C>>( (int)a.degree() );
175        if (a.isConstant() || a.isZERO()) {
176            return roots;
177        }
178        ComplexRing<C> cr = (ComplexRing<C>) a.ring.coFac;
179        // roots of square-free part
180        GenPolynomial<Complex<C>> as = engine.squarefreePart(a);
181        //System.out.println("sqf(a) = " + as);
182        Complex<C> Mb = rootBound(as);
183        C M = Mb.getRe();
184        C M1 = M.sum(M.factory().fromInteger(1)); // asymmetric to origin
185        //System.out.println("M = " + M);
186        if (debug) {
187            logger.info("rootBound = {}", M);
188        }
189        Complex<C>[] corner = (Complex<C>[]) new Complex[4];
190        corner[0] = new Complex<C>(cr, M1.negate(), M); // nw
191        corner[1] = new Complex<C>(cr, M1.negate(), M1.negate()); // sw
192        corner[2] = new Complex<C>(cr, M, M1.negate()); // se
193        corner[3] = new Complex<C>(cr, M, M); // ne
194        Rectangle<C> rect = new Rectangle<C>(corner);
195        List<Rectangle<C>> asr = null;
196        try {
197            asr = complexRoots(rect, as); // must be squarefree
198        } catch (InvalidBoundaryException ex) {
199            //logger.error("invalid boundary {} for as = {}", rect, as);
200            throw new RuntimeException("this should never happen " + ex);
201        }
202        //System.out.println("roots(sqf(a)) = " + asr);
203        // check roots of square-free factors
204        SortedMap<GenPolynomial<Complex<C>>, Long> sa = engine.squarefreeFactors(a); // with multiplicity
205        //System.out.println("squarefree factors = " + sa + ", of a = " + a + ", deg(a) = " + a.degree());
206        for (Map.Entry<GenPolynomial<Complex<C>>, Long> me : sa.entrySet()) {
207            GenPolynomial<Complex<C>> p = me.getKey(); // is squarefree
208            long d = p.degree(0);
209            if (p.isConstant() || p.isZERO()) {
210                continue;
211            }
212            long e = me.getValue();
213            Iterator<Rectangle<C>> rit = asr.iterator();
214            while (rit.hasNext()) {
215                Rectangle<C> rsi = rit.next();
216                long w = 0L;
217                try {
218                    w = complexRootCount(rsi, p); // roots of p not required
219                } catch (InvalidBoundaryException ex) {
220                    //logger.error("invalid boundary {} for p = {}", rsi, p);
221                    throw new RuntimeException("this should never happen " + ex);
222                }
223                if (w == 1) { // since roots of as are separated and must belong to only one p
224                    for (int i = 0; i < e; i++) { // add with multiplicity
225                        roots.add(rsi);
226                    }
227                    rit.remove(); // rsi // removes last object from next()
228                    d--;
229                    if (d <= 0L) {
230                        break; // all roots identified
231                    }
232                }
233            }
234        }
235        //System.out.println("#roots = " + roots.size() + ", roots = " + roots);
236        return roots;
237    }
238
239
240    /**
241     * Complex root refinement of complex polynomial a on rectangle.
242     * @param rect rectangle containing exactly one complex root.
243     * @param a univariate squarefree complex polynomial.
244     * @param len rational length for refinement.
245     * @return refined complex root.
246     */
247    @SuppressWarnings({"cast","unchecked"})
248    public Rectangle<C> complexRootRefinement(Rectangle<C> rect, GenPolynomial<Complex<C>> a, BigRational len)
249        throws InvalidBoundaryException {
250        ComplexRing<C> cr = (ComplexRing<C>) a.ring.coFac;
251        Rectangle<C> root = rect;
252        long w;
253        if (debug) {
254            w = complexRootCount(root, a);
255            if (w != 1) {
256                System.out.println("#root  = " + w + ", root = " + root);
257                System.out.println("deg(a) = " + a.degree() + ", a = " + a);
258                throw new ArithmeticException("no initial isolating rectangle " + rect);
259            }
260        }
261        Complex<C> eps = cr.fromInteger(1);
262        eps = eps.divide(cr.fromInteger(1000)); // 1/1000
263        BigRational length = len.multiply(len);
264        Complex<C> delta = null;
265        boolean work = true;
266        while (work) {
267            try {
268                while (root.rationalLength().compareTo(length) > 0) {
269                    //System.out.println("root = " + root + ", len = " + new BigDecimal(root.rationalLength())); 
270                    if (delta == null) {
271                        delta = root.corners[3].subtract(root.corners[1]);
272                        delta = delta.divide(cr.fromInteger(2));
273                        //System.out.println("delta = " + toDecimal(delta)); 
274                    }
275                    Complex<C> center = root.corners[1].sum(delta);
276                    //System.out.println("refine center = " + toDecimal(center));
277                    if (debug) {
278                        logger.info("new center = {}", center);
279                    }
280
281                    Complex<C>[] cp = (Complex<C>[]) copyOfComplex(root.corners, 4);
282                    // cp[0] fix
283                    cp[1] = new Complex<C>(cr, cp[1].getRe(), center.getIm());
284                    cp[2] = center;
285                    cp[3] = new Complex<C>(cr, center.getRe(), cp[3].getIm());
286                    Rectangle<C> nw = new Rectangle<C>(cp);
287                    w = complexRootCount(nw, a);
288                    if (w == 1) {
289                        root = nw;
290                        delta = null;
291                        continue;
292                    }
293
294                    cp = (Complex<C>[]) copyOfComplex(root.corners, 4);
295                    cp[0] = new Complex<C>(cr, cp[0].getRe(), center.getIm());
296                    // cp[1] fix
297                    cp[2] = new Complex<C>(cr, center.getRe(), cp[2].getIm());
298                    cp[3] = center;
299                    Rectangle<C> sw = new Rectangle<C>(cp);
300                    w = complexRootCount(sw, a);
301                    //System.out.println("#swr = " + w); 
302                    if (w == 1) {
303                        root = sw;
304                        delta = null;
305                        continue;
306                    }
307
308                    cp = (Complex<C>[]) copyOfComplex(root.corners, 4);
309                    cp[0] = center;
310                    cp[1] = new Complex<C>(cr, center.getRe(), cp[1].getIm());
311                    // cp[2] fix
312                    cp[3] = new Complex<C>(cr, cp[3].getRe(), center.getIm());
313                    Rectangle<C> se = new Rectangle<C>(cp);
314                    w = complexRootCount(se, a);
315                    //System.out.println("#ser = " + w); 
316                    if (w == 1) {
317                        root = se;
318                        delta = null;
319                        continue;
320                    }
321
322                    cp = (Complex<C>[]) copyOfComplex(root.corners, 4);
323                    cp[0] = new Complex<C>(cr, center.getRe(), cp[0].getIm());
324                    cp[1] = center;
325                    cp[2] = new Complex<C>(cr, cp[2].getRe(), center.getIm());
326                    // cp[3] fix
327                    Rectangle<C> ne = new Rectangle<C>(cp);
328                    w = complexRootCount(ne, a);
329                    //System.out.println("#ner = " + w); 
330                    if (w == 1) {
331                        root = ne;
332                        delta = null;
333                        continue;
334                    }
335                    if (true) {
336                        w = complexRootCount(root, a);
337                        System.out.println("#root = " + w);
338                        System.out.println("root = " + root);
339                    }
340                    throw new ArithmeticException("no isolating rectangle " + rect);
341                }
342                work = false;
343            } catch (InvalidBoundaryException e) {
344                // repeat with new center
345                delta = delta.sum(delta.multiply(eps)); // distort
346                //System.out.println("new refine delta = " + toDecimal(delta));
347                eps = eps.sum(eps.multiply(cr.getIMAG()));
348            }
349        }
350        return root;
351    }
352
353
354    /**
355     * List of complex roots of complex polynomial.
356     * @param a univariate complex polynomial.
357     * @param len rational length for refinement.
358     * @return list of complex roots to desired precision.
359     */
360    @SuppressWarnings({"cast","unchecked"})
361    public List<Rectangle<C>> complexRoots(GenPolynomial<Complex<C>> a, BigRational len) {
362        ComplexRing<C> cr = (ComplexRing<C>) a.ring.coFac;
363        SortedMap<GenPolynomial<Complex<C>>, Long> sa = engine.squarefreeFactors(a);
364        List<Rectangle<C>> roots = new ArrayList<Rectangle<C>>();
365        for (Map.Entry<GenPolynomial<Complex<C>>, Long> me : sa.entrySet()) {
366            GenPolynomial<Complex<C>> p = me.getKey();
367            Complex<C> Mb = rootBound(p);
368            C M = Mb.getRe();
369            C M1 = M.sum(M.factory().fromInteger(1)); // asymmetric to origin
370            if (debug) {
371                logger.info("rootBound = {}", M);
372            }
373            Complex<C>[] corner = (Complex<C>[]) new Complex[4];
374            corner[0] = new Complex<C>(cr, M1.negate(), M); // nw
375            corner[1] = new Complex<C>(cr, M1.negate(), M1.negate()); // sw
376            corner[2] = new Complex<C>(cr, M, M1.negate()); // se
377            corner[3] = new Complex<C>(cr, M, M); // ne
378            Rectangle<C> rect = new Rectangle<C>(corner);
379            try {
380                List<Rectangle<C>> rs = complexRoots(rect, p);
381                List<Rectangle<C>> rf = new ArrayList<Rectangle<C>>(rs.size());
382                for (Rectangle<C> r : rs) {
383                    Rectangle<C> rr = complexRootRefinement(r, p, len);
384                    rf.add(rr);
385                }
386                long e = me.getValue(); // sa.get(p);
387                for (int i = 0; i < e; i++) { // add with multiplicity
388                    roots.addAll(rf);
389                }
390            } catch (InvalidBoundaryException e) {
391                throw new RuntimeException("this should never happen " + e);
392            }
393        }
394        return roots;
395    }
396
397
398    /**
399     * Invariant rectangle for algebraic number.
400     * @param rect root isolating rectangle for f which contains exactly one
401     *            root.
402     * @param f univariate polynomial, non-zero.
403     * @param g univariate polynomial, gcd(f,g) == 1.
404     * @return v with v a new rectangle contained in iv such that g(w) != 0 for
405     *         w in v.
406     */
407    public abstract Rectangle<C> invariantRectangle(Rectangle<C> rect, GenPolynomial<Complex<C>> f,
408                                                    GenPolynomial<Complex<C>> g) throws InvalidBoundaryException;
409
410
411    /**
412     * Get decimal approximation.
413     * @param a complex number.
414     * @return decimal(a).
415     */
416    public String toDecimal(Complex<C> a) {
417        C r = a.getRe();
418        String s = r.toString();
419        BigRational rs = new BigRational(s);
420        BigDecimal rd = new BigDecimal(rs);
421        C i = a.getIm();
422        s = i.toString();
423        BigRational is = new BigRational(s);
424        BigDecimal id = new BigDecimal(is);
425        //System.out.println("rd = " + rd);
426        //System.out.println("id = " + id);
427        return rd.toString() + " i " + id.toString();
428    }
429
430
431    /**
432     * Approximate complex root.
433     * @param rt root isolating rectangle.
434     * @param f univariate polynomial, non-zero.
435     * @param eps requested interval length.
436     * @return a decimal approximation d such that |d-v| &lt; eps, for f(v) = 0,
437     *         v in rt.
438     */
439    public Complex<BigDecimal> approximateRoot(Rectangle<C> rt, GenPolynomial<Complex<C>> f, BigRational eps)
440        throws NoConvergenceException {
441        if (rt == null) {
442            throw new IllegalArgumentException("null interval not allowed");
443        }
444        Complex<BigDecimal> d = rt.getDecimalCenter();
445        //System.out.println("d  = " + d);
446        if (f == null || f.isZERO() || f.isConstant() || eps == null) {
447            return d;
448        }
449        if (rt.rationalLength().compareTo(eps) < 0) {
450            return d;
451        }
452        ComplexRing<BigDecimal> cr = d.ring;
453        Complex<C> sw = rt.getSW();
454        BigDecimal swr = new BigDecimal(sw.getRe().getRational());
455        BigDecimal swi = new BigDecimal(sw.getIm().getRational());
456        Complex<BigDecimal> ll = new Complex<BigDecimal>(cr, swr, swi);
457        Complex<C> ne = rt.getNE();
458        BigDecimal ner = new BigDecimal(ne.getRe().getRational());
459        BigDecimal nei = new BigDecimal(ne.getIm().getRational());
460        Complex<BigDecimal> ur = new Complex<BigDecimal>(cr, ner, nei);
461
462        BigDecimal e = new BigDecimal(eps.getRational());
463        Complex<BigDecimal> q = new Complex<BigDecimal>(cr, new BigDecimal("0.25"));
464        e = e.multiply(d.norm().getRe()); // relative error
465        //System.out.println("e  = " + e);
466
467        // polynomials with decimal coefficients
468        GenPolynomialRing<Complex<BigDecimal>> dfac = new GenPolynomialRing<Complex<BigDecimal>>(cr, f.ring);
469        GenPolynomial<Complex<BigDecimal>> df = PolyUtil.<C> complexDecimalFromRational(dfac, f);
470        GenPolynomial<Complex<C>> fp = PolyUtil.<Complex<C>> baseDeriviative(f);
471        GenPolynomial<Complex<BigDecimal>> dfp = PolyUtil.<C> complexDecimalFromRational(dfac, fp);
472
473        // Newton Raphson iteration: x_{n+1} = x_n - f(x_n)/f'(x_n)
474        int i = 0;
475        final int MITER = 50;
476        int dir = -1;
477        while (i++ < MITER) {
478            Complex<BigDecimal> fx = PolyUtil.<Complex<BigDecimal>> evaluateMain(cr, df, d); // f(d)
479            //BigDecimal fs = fx.norm().getRe();
480            //System.out.println("fs = " + fs);
481            if (fx.isZERO()) {
482                return d;
483            }
484            Complex<BigDecimal> fpx = PolyUtil.<Complex<BigDecimal>> evaluateMain(cr, dfp, d); // f'(d)
485            if (fpx.isZERO()) {
486                throw new NoConvergenceException("zero deriviative should not happen");
487            }
488            Complex<BigDecimal> x = fx.divide(fpx);
489            Complex<BigDecimal> dx = d.subtract(x);
490            //System.out.println("dx = " + dx);
491            if (d.subtract(dx).norm().getRe().compareTo(e) <= 0) {
492                return dx;
493            }
494            //             if ( false ) { // not useful:
495            //                 Complex<BigDecimal> fxx  = PolyUtil.<Complex<BigDecimal>> evaluateMain(cr, df, dx); // f(dx)
496            //                 //System.out.println("fxx = " + fxx);
497            //                 BigDecimal fsx = fxx.norm().getRe();
498            //                 System.out.println("fsx = " + fsx);
499            //                 while ( fsx.compareTo( fs ) >= 0 ) {
500            //                     System.out.println("trying to increase f(d) ");
501            //                     if ( i++ > MITER ) { // dx > right: dx - right > 0
502            //                         throw new NoConvergenceException("no convergence after " + i + " steps");
503            //                     }
504            //                     x = x.multiply(q); // x * 1/4
505            //                     dx = d.subtract(x);
506            //                     //System.out.println(" x = " + x);
507            //                     System.out.println("dx = " + dx);
508            //                     fxx  = PolyUtil.<Complex<BigDecimal>> evaluateMain(cr, df, dx); // f(dx)
509            //                     //System.out.println("fxx = " + fxx);
510            //                     fsx = fxx.norm().getRe();
511            //                     System.out.println("fsx = " + fsx);
512            //                 }
513            //             }
514            // check interval bounds
515            while (dx.getRe().compareTo(ll.getRe()) < 0 || dx.getIm().compareTo(ll.getIm()) < 0
516                   || dx.getRe().compareTo(ur.getRe()) > 0 || dx.getIm().compareTo(ur.getIm()) > 0) {
517                // dx < ll: dx - ll < 0
518                // dx > ur: dx - ur > 0
519                if (i++ > MITER) { // dx > right: dx - right > 0
520                    throw new NoConvergenceException("no convergence after " + i + " steps");
521                }
522                if (i > MITER / 2 && dir == 0) {
523                    Complex<C> cc = rt.getCenter();
524                    Rectangle<C> nrt = rt.exchangeSE(cc);
525                    Complex<BigDecimal> sd = nrt.getDecimalCenter();
526                    d = sd;
527                    x = cr.getZERO();
528                    logger.info("trying new SE starting point {}", d);
529                    i = 0;
530                    dir = 1;
531                }
532                if (i > MITER / 2 && dir == 1) {
533                    Complex<C> cc = rt.getCenter();
534                    Rectangle<C> nrt = rt.exchangeNW(cc);
535                    Complex<BigDecimal> sd = nrt.getDecimalCenter();
536                    d = sd;
537                    x = cr.getZERO();
538                    logger.info("trying new NW starting point {}", d);
539                    i = 0;
540                    dir = 2;
541                }
542                if (i > MITER / 2 && dir == 2) {
543                    Complex<C> cc = rt.getCenter();
544                    Rectangle<C> nrt = rt.exchangeSW(cc);
545                    Complex<BigDecimal> sd = nrt.getDecimalCenter();
546                    d = sd;
547                    x = cr.getZERO();
548                    logger.info("trying new SW starting point {}", d);
549                    i = 0;
550                    dir = 3;
551                }
552                if (i > MITER / 2 && dir == 3) {
553                    Complex<C> cc = rt.getCenter();
554                    Rectangle<C> nrt = rt.exchangeNE(cc);
555                    Complex<BigDecimal> sd = nrt.getDecimalCenter();
556                    d = sd;
557                    x = cr.getZERO();
558                    logger.info("trying new NE starting point {}", d);
559                    i = 0;
560                    dir = 4;
561                }
562                if (i > MITER / 2 && (dir == -1 || dir == 4 || dir == 5)) {
563                    Complex<C> sr = rt.randomPoint();
564                    BigDecimal srr = new BigDecimal(sr.getRe().getRational());
565                    BigDecimal sri = new BigDecimal(sr.getIm().getRational());
566                    Complex<BigDecimal> sd = new Complex<BigDecimal>(cr, srr, sri);
567                    d = sd;
568                    x = cr.getZERO();
569                    logger.info("trying new random starting point {}", d);
570                    if (dir == -1) {
571                        i = 0;
572                        dir = 0;
573                    } else if (dir == 4) {
574                        i = 0;
575                        dir = 5;
576                    } else {
577                        //i = 0; 
578                        dir = 6; // end
579                    }
580                }
581                x = x.multiply(q); // x * 1/4
582                dx = d.subtract(x);
583                //System.out.println(" x = " + x);
584                //System.out.println("dx = " + dx);
585            }
586            d = dx;
587        }
588        throw new NoConvergenceException("no convergence after " + i + " steps");
589    }
590
591
592    /**
593     * List of decimal approximations of complex roots of complex polynomial.
594     * @param a univariate complex polynomial.
595     * @param eps length for refinement.
596     * @return list of complex decimal root approximations to desired precision.
597     */
598    @SuppressWarnings({"cast","unchecked"})
599    public List<Complex<BigDecimal>> approximateRoots(GenPolynomial<Complex<C>> a, BigRational eps) {
600        ComplexRing<C> cr = (ComplexRing<C>) a.ring.coFac;
601        SortedMap<GenPolynomial<Complex<C>>, Long> sa = engine.squarefreeFactors(a);
602        List<Complex<BigDecimal>> roots = new ArrayList<Complex<BigDecimal>>();
603        for (Map.Entry<GenPolynomial<Complex<C>>, Long> me : sa.entrySet()) {
604            GenPolynomial<Complex<C>> p = me.getKey();
605            List<Complex<BigDecimal>> rf = null;
606            if (p.degree(0) <= 1) {
607                Complex<C> tc = p.trailingBaseCoefficient();
608                tc = tc.negate();
609                BigDecimal rr = new BigDecimal(tc.getRe().getRational());
610                BigDecimal ri = new BigDecimal(tc.getIm().getRational());
611                ComplexRing<BigDecimal> crf = new ComplexRing<BigDecimal>(rr);
612                Complex<BigDecimal> r = new Complex<BigDecimal>(crf, rr, ri);
613                rf = new ArrayList<Complex<BigDecimal>>(1);
614                rf.add(r);
615            } else {
616                Complex<C> Mb = rootBound(p);
617                C M = Mb.getRe();
618                C M1 = M.sum(M.factory().fromInteger(1)); // asymmetric to origin
619                if (debug) {
620                    logger.info("rootBound = {}", M);
621                }
622                Complex<C>[] corner = (Complex<C>[]) new Complex[4];
623                corner[0] = new Complex<C>(cr, M1.negate(), M); // nw
624                corner[1] = new Complex<C>(cr, M1.negate(), M1.negate()); // sw
625                corner[2] = new Complex<C>(cr, M, M1.negate()); // se
626                corner[3] = new Complex<C>(cr, M, M); // ne
627                Rectangle<C> rect = new Rectangle<C>(corner);
628                List<Rectangle<C>> rs = null;
629                try {
630                    rs = complexRoots(rect, p);
631                } catch (InvalidBoundaryException e) {
632                    throw new RuntimeException("this should never happen " + e);
633                }
634                rf = new ArrayList<Complex<BigDecimal>>(rs.size());
635                for (Rectangle<C> r : rs) {
636                    Complex<BigDecimal> rr = null;
637                    while (rr == null) {
638                        try {
639                            rr = approximateRoot(r, p, eps);
640                            rf.add(rr);
641                        } catch (NoConvergenceException e) {
642                            // fall back to exact algorithm
643                            BigRational len = r.rationalLength();
644                            len = len.multiply(new BigRational(1, 1000));
645                            try {
646                                r = complexRootRefinement(r, p, len);
647                                logger.info("fall back rootRefinement = {}", r);
648                                //System.out.println("len = " + len);
649                            } catch (InvalidBoundaryException ee) {
650                                throw new RuntimeException("this should never happen " + ee);
651                            }
652                        }
653                    }
654                }
655            }
656            long e = me.getValue(); // sa.get(p);
657            for (int i = 0; i < e; i++) { // add with multiplicity
658                roots.addAll(rf);
659            }
660        }
661        return roots;
662    }
663
664
665    /**
666     * Copy the specified array.
667     * @param original array.
668     * @param newLength new array length.
669     * @return copy of this.
670     */
671    public Complex[] copyOfComplex(Complex[] original, int newLength) {
672        Complex[] copy = new Complex[newLength];
673        System.arraycopy(original, 0, copy, 0, Math.min(original.length, newLength));
674        return copy;
675    }
676
677
678    /**
679     * Invariant rectangle for algebraic number magnitude.
680     * @param rect root isolating rectangle for f which contains exactly one
681     *            root.
682     * @param f univariate polynomial, non-zero.
683     * @param g univariate polynomial, gcd(f,g) == 1.
684     * @param eps length limit for rectangle length.
685     * @return v with v a new rectangle contained in rect such that |g(a) -
686     *         g(b)| &lt; eps for a, b in v in rect.
687     */
688    public Rectangle<C> invariantMagnitudeRectangle(Rectangle<C> rect, GenPolynomial<Complex<C>> f,
689                                                    GenPolynomial<Complex<C>> g, BigRational eps) throws InvalidBoundaryException {
690        Rectangle<C> v = rect;
691        if (g == null || g.isZERO()) {
692            return v;
693        }
694        if (g.isConstant()) {
695            return v;
696        }
697        if (f == null || f.isZERO() || f.isConstant()) { // ?
698            return v;
699        }
700        GenPolynomial<Complex<C>> gp = PolyUtil.<Complex<C>> baseDeriviative(g);
701        //System.out.println("g  = " + g);
702        //System.out.println("gp = " + gp);
703        BigRational B = magnitudeBound(rect, gp).getRational();
704        //System.out.println("B = " + B + " : " + B.getClass());
705
706        BigRational len = v.rationalLength();
707        BigRational half = new BigRational(1, 2);
708
709        BigRational vlen = v.rationalLength();
710        vlen = vlen.multiply(vlen);
711        //eps = eps.multiply(eps);
712        //System.out.println("v = " + v);
713        //System.out.println("vlen = " + vlen);
714        while (B.multiply(vlen).compareTo(eps) >= 0) { // TODO: test squared
715            len = len.multiply(half);
716            v = complexRootRefinement(v, f, len);
717            //System.out.println("v = " + v);
718            vlen = v.rationalLength();
719            vlen = vlen.multiply(vlen);
720            //System.out.println("vlen = " + vlen);
721        }
722        //System.out.println("vlen = " + vlen);
723        return v;
724    }
725
726
727    /**
728     * Complex algebraic number magnitude.
729     * @param rect root isolating rectangle for f which contains exactly one
730     *            root, with rect such that |g(a) - g(b)| &lt; eps for a, b in
731     *            rect.
732     * @param f univariate polynomial, non-zero.
733     * @param g univariate polynomial, gcd(f,g) == 1.
734     * @return g(rect) .
735     */
736    public Complex<C> complexRectangleMagnitude(Rectangle<C> rect, GenPolynomial<Complex<C>> f,
737                                                GenPolynomial<Complex<C>> g) {
738        if (g.isZERO() || g.isConstant()) {
739            return g.leadingBaseCoefficient();
740        }
741        RingFactory<Complex<C>> cfac = f.ring.coFac;
742        //System.out.println("cfac = " + cfac + " : " + cfac.getClass());
743        Complex<C> c = rect.getCenter();
744        Complex<C> ev = PolyUtil.<Complex<C>> evaluateMain(cfac, g, c);
745        return ev;
746    }
747
748
749    /**
750     * Complex algebraic number magnitude.
751     * @param rect root isolating rectangle for f which contains exactly one
752     *            root, with rect such that |g(a) - g(b)| &lt; eps for a, b in
753     *            rect.
754     * @param f univariate polynomial, non-zero.
755     * @param g univariate polynomial, gcd(f,g) == 1.
756     * @param eps length limit for rectangle length.
757     * @return g(rect) .
758     */
759    public Complex<C> complexMagnitude(Rectangle<C> rect, GenPolynomial<Complex<C>> f,
760                                       GenPolynomial<Complex<C>> g, BigRational eps) throws InvalidBoundaryException {
761        if (g.isZERO() || g.isConstant()) {
762            return g.leadingBaseCoefficient();
763        }
764        Rectangle<C> v = invariantMagnitudeRectangle(rect, f, g, eps);
765        //System.out.println("ref = " + ref);
766        return complexRectangleMagnitude(v, f, g);
767    }
768
769}