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