001/*
002 * $Id: ComplexRootsAbstract.java 4961 2014-10-17 18:59:39Z 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("cast")
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    @SuppressWarnings("cast")
215    public Rectangle<C> complexRootRefinement(Rectangle<C> rect, GenPolynomial<Complex<C>> a, BigRational len)
216                    throws InvalidBoundaryException {
217        ComplexRing<C> cr = (ComplexRing<C>) a.ring.coFac;
218        Rectangle<C> root = rect;
219        long w;
220        if (debug) {
221            w = complexRootCount(root, a);
222            if (w != 1) {
223                System.out.println("#root = " + w);
224                System.out.println("root = " + root);
225                throw new ArithmeticException("no initial isolating rectangle " + rect);
226            }
227        }
228        Complex<C> eps = cr.fromInteger(1);
229        eps = eps.divide(cr.fromInteger(1000)); // 1/1000
230        BigRational length = len.multiply(len);
231        Complex<C> delta = null;
232        boolean work = true;
233        while (work) {
234            try {
235                while (root.rationalLength().compareTo(length) > 0) {
236                    //System.out.println("root = " + root + ", len = " + new BigDecimal(root.rationalLength())); 
237                    if (delta == null) {
238                        delta = root.corners[3].subtract(root.corners[1]);
239                        delta = delta.divide(cr.fromInteger(2));
240                        //System.out.println("delta = " + toDecimal(delta)); 
241                    }
242                    Complex<C> center = root.corners[1].sum(delta);
243                    //System.out.println("refine center = " + toDecimal(center)); 
244                    if (debug) {
245                        logger.info("new center = " + center);
246                    }
247
248                    Complex<C>[] cp = (Complex<C>[]) copyOfComplex(root.corners, 4);
249                    // cp[0] fix
250                    cp[1] = new Complex<C>(cr, cp[1].getRe(), center.getIm());
251                    cp[2] = center;
252                    cp[3] = new Complex<C>(cr, center.getRe(), cp[3].getIm());
253                    Rectangle<C> nw = new Rectangle<C>(cp);
254                    w = complexRootCount(nw, a);
255                    if (w == 1) {
256                        root = nw;
257                        delta = null;
258                        continue;
259                    }
260
261                    cp = (Complex<C>[]) copyOfComplex(root.corners, 4);
262                    cp[0] = new Complex<C>(cr, cp[0].getRe(), center.getIm());
263                    // cp[1] fix
264                    cp[2] = new Complex<C>(cr, center.getRe(), cp[2].getIm());
265                    cp[3] = center;
266                    Rectangle<C> sw = new Rectangle<C>(cp);
267                    w = complexRootCount(sw, a);
268                    //System.out.println("#swr = " + w); 
269                    if (w == 1) {
270                        root = sw;
271                        delta = null;
272                        continue;
273                    }
274
275                    cp = (Complex<C>[]) copyOfComplex(root.corners, 4);
276                    cp[0] = center;
277                    cp[1] = new Complex<C>(cr, center.getRe(), cp[1].getIm());
278                    // cp[2] fix
279                    cp[3] = new Complex<C>(cr, cp[3].getRe(), center.getIm());
280                    Rectangle<C> se = new Rectangle<C>(cp);
281                    w = complexRootCount(se, a);
282                    //System.out.println("#ser = " + w); 
283                    if (w == 1) {
284                        root = se;
285                        delta = null;
286                        continue;
287                    }
288
289                    cp = (Complex<C>[]) copyOfComplex(root.corners, 4);
290                    cp[0] = new Complex<C>(cr, center.getRe(), cp[0].getIm());
291                    cp[1] = center;
292                    cp[2] = new Complex<C>(cr, cp[2].getRe(), center.getIm());
293                    // cp[3] fix
294                    Rectangle<C> ne = new Rectangle<C>(cp);
295                    w = complexRootCount(ne, a);
296                    //System.out.println("#ner = " + w); 
297                    if (w == 1) {
298                        root = ne;
299                        delta = null;
300                        continue;
301                    }
302                    if (true) {
303                        w = complexRootCount(root, a);
304                        System.out.println("#root = " + w);
305                        System.out.println("root = " + root);
306                    }
307                    throw new ArithmeticException("no isolating rectangle " + rect);
308                }
309                work = false;
310            } catch (InvalidBoundaryException e) {
311                // repeat with new center
312                delta = delta.sum(delta.multiply(eps)); // distort
313                //System.out.println("new refine delta = " + toDecimal(delta));
314                eps = eps.sum(eps.multiply(cr.getIMAG()));
315            }
316        }
317        return root;
318    }
319
320
321    /**
322     * List of complex roots of complex polynomial.
323     * @param a univariate complex polynomial.
324     * @param len rational length for refinement.
325     * @return list of complex roots to desired precision.
326     */
327    @SuppressWarnings("cast")
328    public List<Rectangle<C>> complexRoots(GenPolynomial<Complex<C>> a, BigRational len) {
329        ComplexRing<C> cr = (ComplexRing<C>) a.ring.coFac;
330        SortedMap<GenPolynomial<Complex<C>>, Long> sa = engine.squarefreeFactors(a);
331        List<Rectangle<C>> roots = new ArrayList<Rectangle<C>>();
332        for (Map.Entry<GenPolynomial<Complex<C>>, Long> me : sa.entrySet()) {
333            GenPolynomial<Complex<C>> p = me.getKey();
334            Complex<C> Mb = rootBound(p);
335            C M = Mb.getRe();
336            C M1 = M.sum(M.factory().fromInteger(1)); // asymmetric to origin
337            if (debug) {
338                logger.info("rootBound = " + M);
339            }
340            Complex<C>[] corner = (Complex<C>[]) new Complex[4];
341            corner[0] = new Complex<C>(cr, M1.negate(), M); // nw
342            corner[1] = new Complex<C>(cr, M1.negate(), M1.negate()); // sw
343            corner[2] = new Complex<C>(cr, M, M1.negate()); // se
344            corner[3] = new Complex<C>(cr, M, M); // ne
345            Rectangle<C> rect = new Rectangle<C>(corner);
346            try {
347                List<Rectangle<C>> rs = complexRoots(rect, p);
348                List<Rectangle<C>> rf = new ArrayList<Rectangle<C>>(rs.size());
349                for (Rectangle<C> r : rs) {
350                    Rectangle<C> rr = complexRootRefinement(r, p, len);
351                    rf.add(rr);
352                }
353                long e = me.getValue(); // sa.get(p);
354                for (int i = 0; i < e; i++) { // add with multiplicity
355                    roots.addAll(rf);
356                }
357            } catch (InvalidBoundaryException e) {
358                throw new RuntimeException("this should never happen " + e);
359            }
360        }
361        return roots;
362    }
363
364
365    /**
366     * Invariant rectangle for algebraic number.
367     * @param rect root isolating rectangle for f which contains exactly one
368     *            root.
369     * @param f univariate polynomial, non-zero.
370     * @param g univariate polynomial, gcd(f,g) == 1.
371     * @return v with v a new rectangle contained in iv such that g(w) != 0 for
372     *         w in v.
373     */
374    public abstract Rectangle<C> invariantRectangle(Rectangle<C> rect, GenPolynomial<Complex<C>> f,
375                    GenPolynomial<Complex<C>> g) throws InvalidBoundaryException;
376
377
378    /**
379     * Get decimal approximation.
380     * @param a complex number.
381     * @return decimal(a).
382     */
383    public String toDecimal(Complex<C> a) {
384        C r = a.getRe();
385        String s = r.toString();
386        BigRational rs = new BigRational(s);
387        BigDecimal rd = new BigDecimal(rs);
388        C i = a.getIm();
389        s = i.toString();
390        BigRational is = new BigRational(s);
391        BigDecimal id = new BigDecimal(is);
392        //System.out.println("rd = " + rd);
393        //System.out.println("id = " + id);
394        return rd.toString() + " i " + id.toString();
395    }
396
397
398    /**
399     * Approximate complex root.
400     * @param rt root isolating rectangle.
401     * @param f univariate polynomial, non-zero.
402     * @param eps requested interval length.
403     * @return a decimal approximation d such that |d-v| &lt; eps, for f(v) = 0,
404     *         v in rt.
405     */
406    public Complex<BigDecimal> approximateRoot(Rectangle<C> rt, GenPolynomial<Complex<C>> f, C eps)
407                    throws NoConvergenceException {
408        if (rt == null) {
409            throw new IllegalArgumentException("null interval not allowed");
410        }
411        Complex<BigDecimal> d = rt.getDecimalCenter();
412        //System.out.println("d  = " + d);
413        if (f == null || f.isZERO() || f.isConstant() || eps == null) {
414            return d;
415        }
416        if (rt.length().compareTo(eps) < 0) {
417            return d;
418        }
419        ComplexRing<BigDecimal> cr = d.ring;
420        Complex<C> sw = rt.getSW();
421        BigDecimal swr = new BigDecimal(sw.getRe().getRational());
422        BigDecimal swi = new BigDecimal(sw.getIm().getRational());
423        Complex<BigDecimal> ll = new Complex<BigDecimal>(cr, swr, swi);
424        Complex<C> ne = rt.getNE();
425        BigDecimal ner = new BigDecimal(ne.getRe().getRational());
426        BigDecimal nei = new BigDecimal(ne.getIm().getRational());
427        Complex<BigDecimal> ur = new Complex<BigDecimal>(cr, ner, nei);
428
429        BigDecimal e = new BigDecimal(eps.getRational());
430        Complex<BigDecimal> q = new Complex<BigDecimal>(cr, new BigDecimal("0.25"));
431        e = e.multiply(d.norm().getRe()); // relative error
432        //System.out.println("e  = " + e);
433
434        // polynomials with decimal coefficients
435        GenPolynomialRing<Complex<BigDecimal>> dfac = new GenPolynomialRing<Complex<BigDecimal>>(cr, f.ring);
436        GenPolynomial<Complex<BigDecimal>> df = PolyUtil.<C> complexDecimalFromRational(dfac, f);
437        GenPolynomial<Complex<C>> fp = PolyUtil.<Complex<C>> baseDeriviative(f);
438        GenPolynomial<Complex<BigDecimal>> dfp = PolyUtil.<C> complexDecimalFromRational(dfac, fp);
439
440        // Newton Raphson iteration: x_{n+1} = x_n - f(x_n)/f'(x_n)
441        int i = 0;
442        final int MITER = 50;
443        int dir = -1;
444        while (i++ < MITER) {
445            Complex<BigDecimal> fx = PolyUtil.<Complex<BigDecimal>> evaluateMain(cr, df, d); // f(d)
446            //BigDecimal fs = fx.norm().getRe();
447            //System.out.println("fs = " + fs);
448            if (fx.isZERO()) {
449                return d;
450            }
451            Complex<BigDecimal> fpx = PolyUtil.<Complex<BigDecimal>> evaluateMain(cr, dfp, d); // f'(d)
452            if (fpx.isZERO()) {
453                throw new NoConvergenceException("zero deriviative should not happen");
454            }
455            Complex<BigDecimal> x = fx.divide(fpx);
456            Complex<BigDecimal> dx = d.subtract(x);
457            //System.out.println("dx = " + dx);
458            if (d.subtract(dx).norm().getRe().compareTo(e) <= 0) {
459                return dx;
460            }
461            //             if ( false ) { // not useful:
462            //                 Complex<BigDecimal> fxx  = PolyUtil.<Complex<BigDecimal>> evaluateMain(cr, df, dx); // f(dx)
463            //                 //System.out.println("fxx = " + fxx);
464            //                 BigDecimal fsx = fxx.norm().getRe();
465            //                 System.out.println("fsx = " + fsx);
466            //                 while ( fsx.compareTo( fs ) >= 0 ) {
467            //                     System.out.println("trying to increase f(d) ");
468            //                     if ( i++ > MITER ) { // dx > right: dx - right > 0
469            //                         throw new NoConvergenceException("no convergence after " + i + " steps");
470            //                     }
471            //                     x = x.multiply(q); // x * 1/4
472            //                     dx = d.subtract(x);
473            //                     //System.out.println(" x = " + x);
474            //                     System.out.println("dx = " + dx);
475            //                     fxx  = PolyUtil.<Complex<BigDecimal>> evaluateMain(cr, df, dx); // f(dx)
476            //                     //System.out.println("fxx = " + fxx);
477            //                     fsx = fxx.norm().getRe();
478            //                     System.out.println("fsx = " + fsx);
479            //                 }
480            //             }
481            // check interval bounds
482            while (dx.getRe().compareTo(ll.getRe()) < 0 || dx.getIm().compareTo(ll.getIm()) < 0
483                            || dx.getRe().compareTo(ur.getRe()) > 0 || dx.getIm().compareTo(ur.getIm()) > 0) { // dx < ll: dx - ll < 0
484                                                                                                               // dx > ur: dx - ur > 0
485                if (i++ > MITER) { // dx > right: dx - right > 0
486                    throw new NoConvergenceException("no convergence after " + i + " steps");
487                }
488                if (i > MITER / 2 && dir == 0) {
489                    Complex<C> cc = rt.getCenter();
490                    Rectangle<C> nrt = rt.exchangeSE(cc);
491                    Complex<BigDecimal> sd = nrt.getDecimalCenter();
492                    d = sd;
493                    x = cr.getZERO();
494                    logger.info("trying new SE starting point " + d);
495                    i = 0;
496                    dir = 1;
497                }
498                if (i > MITER / 2 && dir == 1) {
499                    Complex<C> cc = rt.getCenter();
500                    Rectangle<C> nrt = rt.exchangeNW(cc);
501                    Complex<BigDecimal> sd = nrt.getDecimalCenter();
502                    d = sd;
503                    x = cr.getZERO();
504                    logger.info("trying new NW starting point " + d);
505                    i = 0;
506                    dir = 2;
507                }
508                if (i > MITER / 2 && dir == 2) {
509                    Complex<C> cc = rt.getCenter();
510                    Rectangle<C> nrt = rt.exchangeSW(cc);
511                    Complex<BigDecimal> sd = nrt.getDecimalCenter();
512                    d = sd;
513                    x = cr.getZERO();
514                    logger.info("trying new SW starting point " + d);
515                    i = 0;
516                    dir = 3;
517                }
518                if (i > MITER / 2 && dir == 3) {
519                    Complex<C> cc = rt.getCenter();
520                    Rectangle<C> nrt = rt.exchangeNE(cc);
521                    Complex<BigDecimal> sd = nrt.getDecimalCenter();
522                    d = sd;
523                    x = cr.getZERO();
524                    logger.info("trying new NE starting point " + d);
525                    i = 0;
526                    dir = 4;
527                }
528                if (i > MITER / 2 && (dir == -1 || dir == 4 || dir == 5)) {
529                    Complex<C> sr = rt.randomPoint();
530                    BigDecimal srr = new BigDecimal(sr.getRe().getRational());
531                    BigDecimal sri = new BigDecimal(sr.getIm().getRational());
532                    Complex<BigDecimal> sd = new Complex<BigDecimal>(cr, srr, sri);
533                    d = sd;
534                    x = cr.getZERO();
535                    logger.info("trying new random starting point " + d);
536                    if (dir == -1) {
537                        i = 0;
538                        dir = 0;
539                    } else if (dir == 4) {
540                        i = 0;
541                        dir = 5;
542                    } else {
543                        //i = 0; 
544                        dir = 6; // end
545                    }
546                }
547                x = x.multiply(q); // x * 1/4
548                dx = d.subtract(x);
549                //System.out.println(" x = " + x);
550                //System.out.println("dx = " + dx);
551            }
552            d = dx;
553        }
554        throw new NoConvergenceException("no convergence after " + i + " steps");
555    }
556
557
558    /**
559     * List of decimal approximations of complex roots of complex polynomial.
560     * @param a univariate complex polynomial.
561     * @param eps length for refinement.
562     * @return list of complex decimal root approximations to desired precision.
563     */
564    @SuppressWarnings("cast")
565    public List<Complex<BigDecimal>> approximateRoots(GenPolynomial<Complex<C>> a, C eps) {
566        ComplexRing<C> cr = (ComplexRing<C>) a.ring.coFac;
567        SortedMap<GenPolynomial<Complex<C>>, Long> sa = engine.squarefreeFactors(a);
568        List<Complex<BigDecimal>> roots = new ArrayList<Complex<BigDecimal>>();
569        for (Map.Entry<GenPolynomial<Complex<C>>, Long> me : sa.entrySet()) {
570            GenPolynomial<Complex<C>> p = me.getKey();
571            List<Complex<BigDecimal>> rf = null;
572            if (p.degree(0) <= 1) {
573                Complex<C> tc = p.trailingBaseCoefficient();
574                tc = tc.negate();
575                BigDecimal rr = new BigDecimal(tc.getRe().getRational());
576                BigDecimal ri = new BigDecimal(tc.getIm().getRational());
577                ComplexRing<BigDecimal> crf = new ComplexRing<BigDecimal>(rr);
578                Complex<BigDecimal> r = new Complex<BigDecimal>(crf, rr, ri);
579                rf = new ArrayList<Complex<BigDecimal>>(1);
580                rf.add(r);
581            } else {
582                Complex<C> Mb = rootBound(p);
583                C M = Mb.getRe();
584                C M1 = M.sum(M.factory().fromInteger(1)); // asymmetric to origin
585                if (debug) {
586                    logger.info("rootBound = " + M);
587                }
588                Complex<C>[] corner = (Complex<C>[]) new Complex[4];
589                corner[0] = new Complex<C>(cr, M1.negate(), M); // nw
590                corner[1] = new Complex<C>(cr, M1.negate(), M1.negate()); // sw
591                corner[2] = new Complex<C>(cr, M, M1.negate()); // se
592                corner[3] = new Complex<C>(cr, M, M); // ne
593                Rectangle<C> rect = new Rectangle<C>(corner);
594                List<Rectangle<C>> rs = null;
595                try {
596                    rs = complexRoots(rect, p);
597                } catch (InvalidBoundaryException e) {
598                    throw new RuntimeException("this should never happen " + e);
599                }
600                rf = new ArrayList<Complex<BigDecimal>>(rs.size());
601                for (Rectangle<C> r : rs) {
602                    Complex<BigDecimal> rr = null;
603                    while (rr == null) {
604                        try {
605                            rr = approximateRoot(r, p, eps);
606                            rf.add(rr);
607                        } catch (NoConvergenceException e) {
608                            // fall back to exact algorithm
609                            BigRational len = r.rationalLength();
610                            len = len.multiply(new BigRational(1, 1000));
611                            try {
612                                r = complexRootRefinement(r, p, len);
613                                logger.info("fall back rootRefinement = " + r);
614                                //System.out.println("len = " + len);
615                            } catch (InvalidBoundaryException ee) {
616                                throw new RuntimeException("this should never happen " + ee);
617                            }
618                        }
619                    }
620                }
621            }
622            long e = me.getValue(); // sa.get(p);
623            for (int i = 0; i < e; i++) { // add with multiplicity
624                roots.addAll(rf);
625            }
626        }
627        return roots;
628    }
629
630
631    /**
632     * Copy the specified array.
633     * @param original array.
634     * @param newLength new array length.
635     * @return copy of this.
636     */
637    public Complex[] copyOfComplex(Complex[] original, int newLength) {
638        Complex[] copy = new Complex[newLength];
639        System.arraycopy(original, 0, copy, 0, Math.min(original.length, newLength));
640        return copy;
641    }
642
643
644    /**
645     * Invariant rectangle for algebraic number magnitude.
646     * @param rect root isolating rectangle for f which contains exactly one
647     *            root.
648     * @param f univariate polynomial, non-zero.
649     * @param g univariate polynomial, gcd(f,g) == 1.
650     * @param eps length limit for rectangle length.
651     * @return v with v a new rectangle contained in rect such that |g(a) -
652     *         g(b)| &lt; eps for a, b in v in rect.
653     */
654    public Rectangle<C> invariantMagnitudeRectangle(Rectangle<C> rect, GenPolynomial<Complex<C>> f,
655                    GenPolynomial<Complex<C>> g, C eps) throws InvalidBoundaryException {
656        Rectangle<C> v = rect;
657        if (g == null || g.isZERO()) {
658            return v;
659        }
660        if (g.isConstant()) {
661            return v;
662        }
663        if (f == null || f.isZERO() || f.isConstant()) { // ?
664            return v;
665        }
666        GenPolynomial<Complex<C>> gp = PolyUtil.<Complex<C>> baseDeriviative(g);
667        //System.out.println("g  = " + g);
668        //System.out.println("gp = " + gp);
669        C B = magnitudeBound(rect, gp);
670        //System.out.println("B = " + B + " : " + B.getClass());
671
672        BigRational len = v.rationalLength();
673        BigRational half = new BigRational(1, 2);
674
675        C vlen = v.length();
676        vlen = vlen.multiply(vlen);
677        //eps = eps.multiply(eps);
678        //System.out.println("v = " + v);
679        //System.out.println("vlen = " + vlen);
680        while (B.multiply(vlen).compareTo(eps) >= 0) { // TODO: test squared
681            len = len.multiply(half);
682            v = complexRootRefinement(v, f, len);
683            //System.out.println("v = " + v);
684            vlen = v.length();
685            vlen = vlen.multiply(vlen);
686            //System.out.println("vlen = " + vlen);
687        }
688        //System.out.println("vlen = " + vlen);
689        return v;
690    }
691
692
693    /**
694     * Complex algebraic number magnitude.
695     * @param rect root isolating rectangle for f which contains exactly one
696     *            root, with rect such that |g(a) - g(b)| &lt; eps for a, b in
697     *            rect.
698     * @param f univariate polynomial, non-zero.
699     * @param g univariate polynomial, gcd(f,g) == 1.
700     * @return g(rect) .
701     */
702    public Complex<C> complexRectangleMagnitude(Rectangle<C> rect, GenPolynomial<Complex<C>> f,
703                    GenPolynomial<Complex<C>> g) {
704        if (g.isZERO() || g.isConstant()) {
705            return g.leadingBaseCoefficient();
706        }
707        RingFactory<Complex<C>> cfac = f.ring.coFac;
708        //System.out.println("cfac = " + cfac + " : " + cfac.getClass());
709        Complex<C> c = rect.getCenter();
710        Complex<C> ev = PolyUtil.<Complex<C>> evaluateMain(cfac, g, c);
711        return ev;
712    }
713
714
715    /**
716     * Complex algebraic number magnitude.
717     * @param rect root isolating rectangle for f which contains exactly one
718     *            root, with rect such that |g(a) - g(b)| &lt; eps for a, b in
719     *            rect.
720     * @param f univariate polynomial, non-zero.
721     * @param g univariate polynomial, gcd(f,g) == 1.
722     * @param eps length limit for rectangle length.
723     * @return g(rect) .
724     */
725    public Complex<C> complexMagnitude(Rectangle<C> rect, GenPolynomial<Complex<C>> f,
726                    GenPolynomial<Complex<C>> g, C eps) throws InvalidBoundaryException {
727        if (g.isZERO() || g.isConstant()) {
728            return g.leadingBaseCoefficient();
729        }
730        Rectangle<C> v = invariantMagnitudeRectangle(rect, f, g, eps);
731        //System.out.println("ref = " + ref);
732        return complexRectangleMagnitude(v, f, g);
733    }
734
735}