View Javadoc
1   /*
2    * Licensed to the Hipparchus project under one or more
3    * contributor license agreements.  See the NOTICE file distributed with
4    * this work for additional information regarding copyright ownership.
5    * The Hipparchus project licenses this file to You under the Apache License, Version 2.0
6    * (the "License"); you may not use this file except in compliance with
7    * the License.  You may obtain a copy of the License at
8    *
9    *      https://www.apache.org/licenses/LICENSE-2.0
10   *
11   * Unless required by applicable law or agreed to in writing, software
12   * distributed under the License is distributed on an "AS IS" BASIS,
13   * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14   * See the License for the specific language governing permissions and
15   * limitations under the License.
16   */
17  package org.hipparchus.complex;
18  
19  import java.util.ArrayList;
20  import java.util.List;
21  
22  import org.hipparchus.CalculusFieldElement;
23  import org.hipparchus.Field;
24  import org.hipparchus.exception.LocalizedCoreFormats;
25  import org.hipparchus.exception.MathIllegalArgumentException;
26  import org.hipparchus.exception.NullArgumentException;
27  import org.hipparchus.util.FastMath;
28  import org.hipparchus.util.FieldSinCos;
29  import org.hipparchus.util.FieldSinhCosh;
30  import org.hipparchus.util.MathArrays;
31  import org.hipparchus.util.MathUtils;
32  import org.hipparchus.util.Precision;
33  
34  /**
35   * Representation of a Complex number, i.e. a number which has both a
36   * real and imaginary part.
37   * <p>
38   * Implementations of arithmetic operations handle {@code NaN} and
39   * infinite values according to the rules for {@link java.lang.Double}, i.e.
40   * {@link #equals} is an equivalence relation for all instances that have
41   * a {@code NaN} in either real or imaginary part, e.g. the following are
42   * considered equal:
43   * <ul>
44   *  <li>{@code 1 + NaNi}</li>
45   *  <li>{@code NaN + i}</li>
46   *  <li>{@code NaN + NaNi}</li>
47   * </ul>
48   * <p>
49   * Note that this contradicts the IEEE-754 standard for floating
50   * point numbers (according to which the test {@code x == x} must fail if
51   * {@code x} is {@code NaN}). The method
52   * {@link org.hipparchus.util.Precision#equals(double,double,int)
53   * equals for primitive double} in {@link org.hipparchus.util.Precision}
54   * conforms with IEEE-754 while this class conforms with the standard behavior
55   * for Java object types.
56   * @param <T> the type of the field elements
57   * @since 2.0
58   */
59  public class FieldComplex<T extends CalculusFieldElement<T>> implements CalculusFieldElement<FieldComplex<T>>  {
60  
61      /** A real number representing log(10). */
62      private static final double LOG10 = 2.302585092994045684;
63  
64      /** The imaginary part. */
65      private final T imaginary;
66  
67      /** The real part. */
68      private final T real;
69  
70      /** Record whether this complex number is equal to NaN. */
71      private final transient boolean isNaN;
72  
73      /** Record whether this complex number is infinite. */
74      private final transient boolean isInfinite;
75  
76      /**
77       * Create a complex number given only the real part.
78       *
79       * @param real Real part.
80       */
81      public FieldComplex(T real) {
82          this(real, real.getField().getZero());
83      }
84  
85      /**
86       * Create a complex number given the real and imaginary parts.
87       *
88       * @param real Real part.
89       * @param imaginary Imaginary part.
90       */
91      public FieldComplex(T real, T imaginary) {
92          this.real = real;
93          this.imaginary = imaginary;
94  
95          isNaN = real.isNaN() || imaginary.isNaN();
96          isInfinite = !isNaN &&
97              (real.isInfinite() || imaginary.isInfinite());
98      }
99  
100     /** Get the square root of -1.
101      * @param field field the complex components belong to
102      * @return number representing "0.0 + 1.0i"
103      * @param <T> the type of the field elements
104      */
105     public static <T extends CalculusFieldElement<T>> FieldComplex<T> getI(final Field<T> field) {
106         return new FieldComplex<>(field.getZero(), field.getOne());
107     }
108 
109     /** Get the square root of -1.
110      * @param field field the complex components belong to
111      * @return number representing "0.0 _ 1.0i"
112      * @param <T> the type of the field elements
113      */
114     public static <T extends CalculusFieldElement<T>> FieldComplex<T> getMinusI(final Field<T> field) {
115         return new FieldComplex<>(field.getZero(), field.getOne().negate());
116     }
117 
118     /** Get a complex number representing "NaN + NaNi".
119      * @param field field the complex components belong to
120      * @return complex number representing "NaN + NaNi"
121      * @param <T> the type of the field elements
122      */
123     public static <T extends CalculusFieldElement<T>> FieldComplex<T> getNaN(final Field<T> field) {
124         return new FieldComplex<>(field.getZero().add(Double.NaN), field.getZero().add(Double.NaN));
125     }
126 
127     /** Get a complex number representing "+INF + INFi".
128      * @param field field the complex components belong to
129      * @return complex number representing "+INF + INFi"
130      * @param <T> the type of the field elements
131      */
132     public static <T extends CalculusFieldElement<T>> FieldComplex<T> getInf(final Field<T> field) {
133         return new FieldComplex<>(field.getZero().add(Double.POSITIVE_INFINITY), field.getZero().add(Double.POSITIVE_INFINITY));
134     }
135 
136     /** Get a complex number representing "1.0 + 0.0i".
137      * @param field field the complex components belong to
138      * @return complex number representing "1.0 + 0.0i"
139      * @param <T> the type of the field elements
140      */
141     public static <T extends CalculusFieldElement<T>> FieldComplex<T> getOne(final Field<T> field) {
142         return new FieldComplex<>(field.getOne(), field.getZero());
143     }
144 
145     /** Get a complex number representing "-1.0 + 0.0i".
146      * @param field field the complex components belong to
147      * @return complex number representing "-1.0 + 0.0i"
148      * @param <T> the type of the field elements
149      */
150     public static <T extends CalculusFieldElement<T>> FieldComplex<T> getMinusOne(final Field<T> field) {
151         return new FieldComplex<>(field.getOne().negate(), field.getZero());
152     }
153 
154     /** Get a complex number representing "0.0 + 0.0i".
155      * @param field field the complex components belong to
156      * @return complex number representing "0.0 + 0.0i
157      * @param <T> the type of the field elements
158      */
159     public static <T extends CalculusFieldElement<T>> FieldComplex<T> getZero(final Field<T> field) {
160         return new FieldComplex<>(field.getZero(), field.getZero());
161     }
162 
163     /** Get a complex number representing "π + 0.0i".
164      * @param field field the complex components belong to
165      * @return complex number representing "π + 0.0i
166      * @param <T> the type of the field elements
167      */
168     public static <T extends CalculusFieldElement<T>> FieldComplex<T> getPi(final Field<T> field) {
169         return new FieldComplex<>(field.getZero().getPi(), field.getZero());
170     }
171 
172     /**
173      * Return the absolute value of this complex number.
174      * Returns {@code NaN} if either real or imaginary part is {@code NaN}
175      * and {@code Double.POSITIVE_INFINITY} if neither part is {@code NaN},
176      * but at least one part is infinite.
177      *
178      * @return the absolute value.
179      */
180     @Override
181     public FieldComplex<T> abs() {
182         // we check NaN here because FastMath.hypot checks it after infinity
183         return isNaN ? getNaN(getPartsField()) : createComplex(FastMath.hypot(real, imaginary), getPartsField().getZero());
184     }
185 
186     /**
187      * Returns a {@code Complex} whose value is
188      * {@code (this + addend)}.
189      * Uses the definitional formula
190      * <p>
191      *   {@code (a + bi) + (c + di) = (a+c) + (b+d)i}
192      * </p>
193      * If either {@code this} or {@code addend} has a {@code NaN} value in
194      * either part, {@link #getNaN(Field)} is returned; otherwise {@code Infinite}
195      * and {@code NaN} values are returned in the parts of the result
196      * according to the rules for {@link java.lang.Double} arithmetic.
197      *
198      * @param  addend Value to be added to this {@code Complex}.
199      * @return {@code this + addend}.
200      * @throws NullArgumentException if {@code addend} is {@code null}.
201      */
202     @Override
203     public FieldComplex<T> add(FieldComplex<T> addend) throws NullArgumentException {
204         MathUtils.checkNotNull(addend);
205         if (isNaN || addend.isNaN) {
206             return getNaN(getPartsField());
207         }
208 
209         return createComplex(real.add(addend.getRealPart()),
210                              imaginary.add(addend.getImaginaryPart()));
211     }
212 
213     /**
214      * Returns a {@code Complex} whose value is {@code (this + addend)},
215      * with {@code addend} interpreted as a real number.
216      *
217      * @param addend Value to be added to this {@code Complex}.
218      * @return {@code this + addend}.
219      * @see #add(FieldComplex)
220      */
221     public FieldComplex<T> add(T addend) {
222         if (isNaN || addend.isNaN()) {
223             return getNaN(getPartsField());
224         }
225 
226         return createComplex(real.add(addend), imaginary);
227     }
228 
229     /**
230      * Returns a {@code Complex} whose value is {@code (this + addend)},
231      * with {@code addend} interpreted as a real number.
232      *
233      * @param addend Value to be added to this {@code Complex}.
234      * @return {@code this + addend}.
235      * @see #add(FieldComplex)
236      */
237     @Override
238     public FieldComplex<T> add(double addend) {
239         if (isNaN || Double.isNaN(addend)) {
240             return getNaN(getPartsField());
241         }
242 
243         return createComplex(real.add(addend), imaginary);
244     }
245 
246     /**
247      * Returns the conjugate of this complex number.
248      * The conjugate of {@code a + bi} is {@code a - bi}.
249      * <p>
250      * {@link #getNaN(Field)} is returned if either the real or imaginary
251      * part of this Complex number equals {@code Double.NaN}.
252      * </p><p>
253      * If the imaginary part is infinite, and the real part is not
254      * {@code NaN}, the returned value has infinite imaginary part
255      * of the opposite sign, e.g. the conjugate of
256      * {@code 1 + POSITIVE_INFINITY i} is {@code 1 - NEGATIVE_INFINITY i}.
257      * </p>
258      * @return the conjugate of this Complex object.
259      */
260     public FieldComplex<T> conjugate() {
261         if (isNaN) {
262             return getNaN(getPartsField());
263         }
264 
265         return createComplex(real, imaginary.negate());
266     }
267 
268     /**
269      * Returns a {@code Complex} whose value is
270      * {@code (this / divisor)}.
271      * Implements the definitional formula
272      * <pre>
273      *  <code>
274      *    a + bi          ac + bd + (bc - ad)i
275      *    ----------- = -------------------------
276      *    c + di         c<sup>2</sup> + d<sup>2</sup>
277      *  </code>
278      * </pre>
279      * but uses
280      * <a href="http://doi.acm.org/10.1145/1039813.1039814">
281      * prescaling of operands</a> to limit the effects of overflows and
282      * underflows in the computation.
283      * <p>
284      * {@code Infinite} and {@code NaN} values are handled according to the
285      * following rules, applied in the order presented:
286      * <ul>
287      *  <li>If either {@code this} or {@code divisor} has a {@code NaN} value
288      *   in either part, {@link #getNaN(Field)} is returned.
289      *  </li>
290      *  <li>If {@code divisor} equals {@link #getZero(Field)}, {@link #getNaN(Field)} is returned.
291      *  </li>
292      *  <li>If {@code this} and {@code divisor} are both infinite,
293      *   {@link #getNaN(Field)} is returned.
294      *  </li>
295      *  <li>If {@code this} is finite (i.e., has no {@code Infinite} or
296      *   {@code NaN} parts) and {@code divisor} is infinite (one or both parts
297      *   infinite), {@link #getZero(Field)} is returned.
298      *  </li>
299      *  <li>If {@code this} is infinite and {@code divisor} is finite,
300      *   {@code NaN} values are returned in the parts of the result if the
301      *   {@link java.lang.Double} rules applied to the definitional formula
302      *   force {@code NaN} results.
303      *  </li>
304      * </ul>
305      *
306      * @param divisor Value by which this {@code Complex} is to be divided.
307      * @return {@code this / divisor}.
308      * @throws NullArgumentException if {@code divisor} is {@code null}.
309      */
310     @Override
311     public FieldComplex<T> divide(FieldComplex<T> divisor)
312         throws NullArgumentException {
313         MathUtils.checkNotNull(divisor);
314         if (isNaN || divisor.isNaN) {
315             return getNaN(getPartsField());
316         }
317 
318         final T c = divisor.getRealPart();
319         final T d = divisor.getImaginaryPart();
320         if (c.isZero() && d.isZero()) {
321             return getNaN(getPartsField());
322         }
323 
324         if (divisor.isInfinite() && !isInfinite()) {
325             return getZero(getPartsField());
326         }
327 
328         if (FastMath.abs(c).getReal() < FastMath.abs(d).getReal()) {
329             T q = c.divide(d);
330             T invDen = c.multiply(q).add(d).reciprocal();
331             return createComplex(real.multiply(q).add(imaginary).multiply(invDen),
332                                  imaginary.multiply(q).subtract(real).multiply(invDen));
333         } else {
334             T q = d.divide(c);
335             T invDen = d.multiply(q).add(c).reciprocal();
336             return createComplex(imaginary.multiply(q).add(real).multiply(invDen),
337                                  imaginary.subtract(real.multiply(q)).multiply(invDen));
338         }
339     }
340 
341     /**
342      * Returns a {@code Complex} whose value is {@code (this / divisor)},
343      * with {@code divisor} interpreted as a real number.
344      *
345      * @param  divisor Value by which this {@code Complex} is to be divided.
346      * @return {@code this / divisor}.
347      * @see #divide(FieldComplex)
348      */
349     public FieldComplex<T> divide(T divisor) {
350         if (isNaN || divisor.isNaN()) {
351             return getNaN(getPartsField());
352         }
353         if (divisor.isZero()) {
354             return getNaN(getPartsField());
355         }
356         if (divisor.isInfinite()) {
357             return !isInfinite() ? getZero(getPartsField()) : getNaN(getPartsField());
358         }
359         return createComplex(real.divide(divisor), imaginary.divide(divisor));
360     }
361 
362     /**
363      * Returns a {@code Complex} whose value is {@code (this / divisor)},
364      * with {@code divisor} interpreted as a real number.
365      *
366      * @param  divisor Value by which this {@code Complex} is to be divided.
367      * @return {@code this / divisor}.
368      * @see #divide(FieldComplex)
369      */
370     @Override
371     public FieldComplex<T> divide(double divisor) {
372         if (isNaN || Double.isNaN(divisor)) {
373             return getNaN(getPartsField());
374         }
375         if (divisor == 0.0) {
376             return getNaN(getPartsField());
377         }
378         if (Double.isInfinite(divisor)) {
379             return !isInfinite() ? getZero(getPartsField()) : getNaN(getPartsField());
380         }
381         return createComplex(real.divide(divisor), imaginary.divide(divisor));
382     }
383 
384     /** {@inheritDoc} */
385     @Override
386     public FieldComplex<T> reciprocal() {
387         if (isNaN) {
388             return getNaN(getPartsField());
389         }
390 
391         if (real.isZero() && imaginary.isZero()) {
392             return getInf(getPartsField());
393         }
394 
395         if (isInfinite) {
396             return getZero(getPartsField());
397         }
398 
399         if (FastMath.abs(real).getReal() < FastMath.abs(imaginary).getReal()) {
400             T q = real.divide(imaginary);
401             T scale = real.multiply(q).add(imaginary).reciprocal();
402             return createComplex(scale.multiply(q), scale.negate());
403         } else {
404             T q = imaginary.divide(real);
405             T scale = imaginary.multiply(q).add(real).reciprocal();
406             return createComplex(scale, scale.negate().multiply(q));
407         }
408     }
409 
410     /**
411      * Test for equality with another object.
412      * If both the real and imaginary parts of two complex numbers
413      * are exactly the same, and neither is {@code Double.NaN}, the two
414      * Complex objects are considered to be equal.
415      * The behavior is the same as for JDK's {@link Double#equals(Object)
416      * Double}:
417      * <ul>
418      *  <li>All {@code NaN} values are considered to be equal,
419      *   i.e, if either (or both) real and imaginary parts of the complex
420      *   number are equal to {@code Double.NaN}, the complex number is equal
421      *   to {@code NaN}.
422      *  </li>
423      *  <li>
424      *   Instances constructed with different representations of zero (i.e.
425      *   either "0" or "-0") are <em>not</em> considered to be equal.
426      *  </li>
427      * </ul>
428      *
429      * @param other Object to test for equality with this instance.
430      * @return {@code true} if the objects are equal, {@code false} if object
431      * is {@code null}, not an instance of {@code Complex}, or not equal to
432      * this instance.
433      */
434     @Override
435     public boolean equals(Object other) {
436         if (this == other) {
437             return true;
438         }
439         if (other instanceof FieldComplex){
440             @SuppressWarnings("unchecked")
441             FieldComplex<T> c = (FieldComplex<T>) other;
442             if (c.isNaN) {
443                 return isNaN;
444             } else {
445                 return real.equals(c.real) && imaginary.equals(c.imaginary);
446             }
447         }
448         return false;
449     }
450 
451     /**
452      * Test for the floating-point equality between Complex objects.
453      * It returns {@code true} if both arguments are equal or within the
454      * range of allowed error (inclusive).
455      *
456      * @param x First value (cannot be {@code null}).
457      * @param y Second value (cannot be {@code null}).
458      * @param maxUlps {@code (maxUlps - 1)} is the number of floating point
459      * values between the real (resp. imaginary) parts of {@code x} and
460      * {@code y}.
461      * @param <T> the type of the field elements
462      * @return {@code true} if there are fewer than {@code maxUlps} floating
463      * point values between the real (resp. imaginary) parts of {@code x}
464      * and {@code y}.
465      *
466      * @see Precision#equals(double,double,int)
467      */
468     public static <T extends CalculusFieldElement<T>>boolean equals(FieldComplex<T> x, FieldComplex<T> y, int maxUlps) {
469         return Precision.equals(x.real.getReal(), y.real.getReal(), maxUlps) &&
470                Precision.equals(x.imaginary.getReal(), y.imaginary.getReal(), maxUlps);
471     }
472 
473     /**
474      * Returns {@code true} iff the values are equal as defined by
475      * {@link #equals(FieldComplex,FieldComplex,int) equals(x, y, 1)}.
476      *
477      * @param x First value (cannot be {@code null}).
478      * @param y Second value (cannot be {@code null}).
479      * @param <T> the type of the field elements
480      * @return {@code true} if the values are equal.
481      */
482     public static <T extends CalculusFieldElement<T>>boolean equals(FieldComplex<T> x, FieldComplex<T> y) {
483         return equals(x, y, 1);
484     }
485 
486     /**
487      * Returns {@code true} if, both for the real part and for the imaginary
488      * part, there is no T value strictly between the arguments or the
489      * difference between them is within the range of allowed error
490      * (inclusive).  Returns {@code false} if either of the arguments is NaN.
491      *
492      * @param x First value (cannot be {@code null}).
493      * @param y Second value (cannot be {@code null}).
494      * @param eps Amount of allowed absolute error.
495      * @param <T> the type of the field elements
496      * @return {@code true} if the values are two adjacent floating point
497      * numbers or they are within range of each other.
498      *
499      * @see Precision#equals(double,double,double)
500      */
501     public static <T extends CalculusFieldElement<T>>boolean equals(FieldComplex<T> x, FieldComplex<T> y,
502                                                                     double eps) {
503         return Precision.equals(x.real.getReal(), y.real.getReal(), eps) &&
504                Precision.equals(x.imaginary.getReal(), y.imaginary.getReal(), eps);
505     }
506 
507     /**
508      * Returns {@code true} if, both for the real part and for the imaginary
509      * part, there is no T value strictly between the arguments or the
510      * relative difference between them is smaller or equal to the given
511      * tolerance. Returns {@code false} if either of the arguments is NaN.
512      *
513      * @param x First value (cannot be {@code null}).
514      * @param y Second value (cannot be {@code null}).
515      * @param eps Amount of allowed relative error.
516      * @param <T> the type of the field elements
517      * @return {@code true} if the values are two adjacent floating point
518      * numbers or they are within range of each other.
519      *
520      * @see Precision#equalsWithRelativeTolerance(double,double,double)
521      */
522     public static <T extends CalculusFieldElement<T>>boolean equalsWithRelativeTolerance(FieldComplex<T> x,
523                                                                                          FieldComplex<T> y,
524                                                                                          double eps) {
525         return Precision.equalsWithRelativeTolerance(x.real.getReal(), y.real.getReal(), eps) &&
526                Precision.equalsWithRelativeTolerance(x.imaginary.getReal(), y.imaginary.getReal(), eps);
527     }
528 
529     /**
530      * Get a hashCode for the complex number.
531      * Any {@code Double.NaN} value in real or imaginary part produces
532      * the same hash code {@code 7}.
533      *
534      * @return a hash code value for this object.
535      */
536     @Override
537     public int hashCode() {
538         if (isNaN) {
539             return 7;
540         }
541         return 37 * (17 * imaginary.hashCode() + real.hashCode());
542     }
543 
544     /** {@inheritDoc}
545      * <p>
546      * This implementation considers +0.0 and -0.0 to be equal for both
547      * real and imaginary components.
548      * </p>
549      */
550     @Override
551     public boolean isZero() {
552         return real.isZero() && imaginary.isZero();
553     }
554 
555     /**
556      * Access the imaginary part.
557      *
558      * @return the imaginary part.
559      */
560     public T getImaginary() {
561         return imaginary;
562     }
563 
564     /**
565      * Access the imaginary part.
566      *
567      * @return the imaginary part.
568      */
569     public T getImaginaryPart() {
570         return imaginary;
571     }
572 
573     /**
574      * Access the real part.
575      *
576      * @return the real part.
577      */
578     @Override
579     public double getReal() {
580         return real.getReal();
581     }
582 
583     /** {@inheritDoc} */
584     @Override
585     public FieldComplex<T> getAddendum() {
586         return new FieldComplex<>(real.getAddendum(), imaginary);
587     }
588 
589     /**
590      * Access the real part.
591      *
592      * @return the real part.
593      */
594     public T getRealPart() {
595         return real;
596     }
597 
598     /**
599      * Checks whether either or both parts of this complex number is
600      * {@code NaN}.
601      *
602      * @return true if either or both parts of this complex number is
603      * {@code NaN}; false otherwise.
604      */
605     @Override
606     public boolean isNaN() {
607         return isNaN;
608     }
609 
610     /** Check whether the instance is real (i.e. imaginary part is zero).
611      * @return true if imaginary part is zero
612       */
613     public boolean isReal() {
614         return imaginary.isZero();
615     }
616 
617     /** Check whether the instance is an integer (i.e. imaginary part is zero and real part has no fractional part).
618      * @return true if imaginary part is zero and real part has no fractional part
619      */
620     public boolean isMathematicalInteger() {
621         return isReal() && Precision.isMathematicalInteger(real.getReal());
622     }
623 
624     /**
625      * Checks whether either the real or imaginary part of this complex number
626      * takes an infinite value (either {@code Double.POSITIVE_INFINITY} or
627      * {@code Double.NEGATIVE_INFINITY}) and neither part
628      * is {@code NaN}.
629      *
630      * @return true if one or both parts of this complex number are infinite
631      * and neither part is {@code NaN}.
632      */
633     @Override
634     public boolean isInfinite() {
635         return isInfinite;
636     }
637 
638     /**
639      * Returns a {@code Complex} whose value is {@code this * factor}.
640      * Implements preliminary checks for {@code NaN} and infinity followed by
641      * the definitional formula:
642      * <p>
643      *   {@code (a + bi)(c + di) = (ac - bd) + (ad + bc)i}
644      * </p>
645      * Returns {@link #getNaN(Field)} if either {@code this} or {@code factor} has one or
646      * more {@code NaN} parts.
647      * <p>
648      * Returns {@link #getInf(Field)} if neither {@code this} nor {@code factor} has one
649      * or more {@code NaN} parts and if either {@code this} or {@code factor}
650      * has one or more infinite parts (same result is returned regardless of
651      * the sign of the components).
652      * </p><p>
653      * Returns finite values in components of the result per the definitional
654      * formula in all remaining cases.</p>
655      *
656      * @param  factor value to be multiplied by this {@code Complex}.
657      * @return {@code this * factor}.
658      * @throws NullArgumentException if {@code factor} is {@code null}.
659      */
660     @Override
661     public FieldComplex<T> multiply(FieldComplex<T> factor)
662         throws NullArgumentException {
663         MathUtils.checkNotNull(factor);
664         if (isNaN || factor.isNaN) {
665             return getNaN(getPartsField());
666         }
667         if (real.isInfinite() ||
668             imaginary.isInfinite() ||
669             factor.real.isInfinite() ||
670             factor.imaginary.isInfinite()) {
671             // we don't use isInfinite() to avoid testing for NaN again
672             return getInf(getPartsField());
673         }
674         return createComplex(real.linearCombination(real, factor.real, imaginary.negate(), factor.imaginary),
675                              real.linearCombination(real, factor.imaginary, imaginary, factor.real));
676     }
677 
678     /**
679      * Returns a {@code Complex} whose value is {@code this * factor}, with {@code factor}
680      * interpreted as a integer number.
681      *
682      * @param  factor value to be multiplied by this {@code Complex}.
683      * @return {@code this * factor}.
684      * @see #multiply(FieldComplex)
685      */
686     @Override
687     public FieldComplex<T> multiply(final int factor) {
688         if (isNaN) {
689             return getNaN(getPartsField());
690         }
691         if (real.isInfinite() || imaginary.isInfinite()) {
692             return getInf(getPartsField());
693         }
694         return createComplex(real.multiply(factor), imaginary.multiply(factor));
695     }
696 
697     /**
698      * Returns a {@code Complex} whose value is {@code this * factor}, with {@code factor}
699      * interpreted as a real number.
700      *
701      * @param  factor value to be multiplied by this {@code Complex}.
702      * @return {@code this * factor}.
703      * @see #multiply(FieldComplex)
704      */
705     @Override
706     public FieldComplex<T> multiply(double factor) {
707         if (isNaN || Double.isNaN(factor)) {
708             return getNaN(getPartsField());
709         }
710         if (real.isInfinite() ||
711             imaginary.isInfinite() ||
712             Double.isInfinite(factor)) {
713             // we don't use isInfinite() to avoid testing for NaN again
714             return getInf(getPartsField());
715         }
716         return createComplex(real.multiply(factor), imaginary.multiply(factor));
717     }
718 
719     /**
720      * Returns a {@code Complex} whose value is {@code this * factor}, with {@code factor}
721      * interpreted as a real number.
722      *
723      * @param  factor value to be multiplied by this {@code Complex}.
724      * @return {@code this * factor}.
725      * @see #multiply(FieldComplex)
726      */
727     public FieldComplex<T> multiply(T factor) {
728         if (isNaN || factor.isNaN()) {
729             return getNaN(getPartsField());
730         }
731         if (real.isInfinite() ||
732             imaginary.isInfinite() ||
733             factor.isInfinite()) {
734             // we don't use isInfinite() to avoid testing for NaN again
735             return getInf(getPartsField());
736         }
737         return createComplex(real.multiply(factor), imaginary.multiply(factor));
738     }
739 
740     /** Compute this * i.
741      * @return this * i
742      * @since 2.0
743      */
744     public FieldComplex<T> multiplyPlusI() {
745         return createComplex(imaginary.negate(), real);
746     }
747 
748     /** Compute this *- -i.
749      * @return this * i
750      * @since 2.0
751      */
752     public FieldComplex<T> multiplyMinusI() {
753         return createComplex(imaginary, real.negate());
754     }
755 
756     @Override
757     public FieldComplex<T> square() {
758         return multiply(this);
759     }
760 
761     /**
762      * Returns a {@code Complex} whose value is {@code (-this)}.
763      * Returns {@code NaN} if either real or imaginary
764      * part of this Complex number is {@code Double.NaN}.
765      *
766      * @return {@code -this}.
767      */
768     @Override
769     public FieldComplex<T> negate() {
770         if (isNaN) {
771             return getNaN(getPartsField());
772         }
773 
774         return createComplex(real.negate(), imaginary.negate());
775     }
776 
777     /**
778      * Returns a {@code Complex} whose value is
779      * {@code (this - subtrahend)}.
780      * Uses the definitional formula
781      * <p>
782      *  {@code (a + bi) - (c + di) = (a-c) + (b-d)i}
783      * </p>
784      * If either {@code this} or {@code subtrahend} has a {@code NaN]} value in either part,
785      * {@link #getNaN(Field)} is returned; otherwise infinite and {@code NaN} values are
786      * returned in the parts of the result according to the rules for
787      * {@link java.lang.Double} arithmetic.
788      *
789      * @param  subtrahend value to be subtracted from this {@code Complex}.
790      * @return {@code this - subtrahend}.
791      * @throws NullArgumentException if {@code subtrahend} is {@code null}.
792      */
793     @Override
794     public FieldComplex<T> subtract(FieldComplex<T> subtrahend)
795         throws NullArgumentException {
796         MathUtils.checkNotNull(subtrahend);
797         if (isNaN || subtrahend.isNaN) {
798             return getNaN(getPartsField());
799         }
800 
801         return createComplex(real.subtract(subtrahend.getRealPart()),
802                              imaginary.subtract(subtrahend.getImaginaryPart()));
803     }
804 
805     /**
806      * Returns a {@code Complex} whose value is
807      * {@code (this - subtrahend)}.
808      *
809      * @param  subtrahend value to be subtracted from this {@code Complex}.
810      * @return {@code this - subtrahend}.
811      * @see #subtract(FieldComplex)
812      */
813     @Override
814     public FieldComplex<T> subtract(double subtrahend) {
815         if (isNaN || Double.isNaN(subtrahend)) {
816             return getNaN(getPartsField());
817         }
818         return createComplex(real.subtract(subtrahend), imaginary);
819     }
820 
821     /**
822      * Returns a {@code Complex} whose value is
823      * {@code (this - subtrahend)}.
824      *
825      * @param  subtrahend value to be subtracted from this {@code Complex}.
826      * @return {@code this - subtrahend}.
827      * @see #subtract(FieldComplex)
828      */
829     public FieldComplex<T> subtract(T subtrahend) {
830         if (isNaN || subtrahend.isNaN()) {
831             return getNaN(getPartsField());
832         }
833         return createComplex(real.subtract(subtrahend), imaginary);
834     }
835 
836     /**
837      * Compute the
838      * <a href="http://mathworld.wolfram.com/InverseCosine.html" TARGET="_top">
839      * inverse cosine</a> of this complex number.
840      * Implements the formula:
841      * <p>
842      *  {@code acos(z) = -i (log(z + i (sqrt(1 - z<sup>2</sup>))))}
843      * </p>
844      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
845      * input argument is {@code NaN} or infinite.
846      *
847      * @return the inverse cosine of this complex number.
848      */
849     @Override
850     public FieldComplex<T> acos() {
851         if (isNaN) {
852             return getNaN(getPartsField());
853         }
854 
855         return this.add(this.sqrt1z().multiplyPlusI()).log().multiplyMinusI();
856     }
857 
858     /**
859      * Compute the
860      * <a href="http://mathworld.wolfram.com/InverseSine.html" TARGET="_top">
861      * inverse sine</a> of this complex number.
862      * Implements the formula:
863      * <p>
864      *  {@code asin(z) = -i (log(sqrt(1 - z<sup>2</sup>) + iz))}
865      * </p><p>
866      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
867      * input argument is {@code NaN} or infinite.</p>
868      *
869      * @return the inverse sine of this complex number.
870      */
871     @Override
872     public FieldComplex<T> asin() {
873         if (isNaN) {
874             return getNaN(getPartsField());
875         }
876 
877         return sqrt1z().add(this.multiplyPlusI()).log().multiplyMinusI();
878     }
879 
880     /**
881      * Compute the
882      * <a href="http://mathworld.wolfram.com/InverseTangent.html" TARGET="_top">
883      * inverse tangent</a> of this complex number.
884      * Implements the formula:
885      * <p>
886      * {@code atan(z) = (i/2) log((1 - iz)/(1 + iz))}
887      * </p><p>
888      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
889      * input argument is {@code NaN} or infinite.</p>
890      *
891      * @return the inverse tangent of this complex number
892      */
893     @Override
894     public FieldComplex<T> atan() {
895         if (isNaN) {
896             return getNaN(getPartsField());
897         }
898 
899         final T one = getPartsField().getOne();
900         if (real.isZero()) {
901 
902             // singularity at ±i
903             if (imaginary.multiply(imaginary).subtract(one).isZero()) {
904                 return getNaN(getPartsField());
905             }
906 
907             // branch cut on imaginary axis
908             final T zero = getPartsField().getZero();
909             final FieldComplex<T> tmp = createComplex(one.add(imaginary).divide(one.subtract(imaginary)), zero).
910                                         log().multiplyPlusI().multiply(0.5);
911             return createComplex(FastMath.copySign(tmp.real, real), tmp.imaginary);
912 
913         } else if (imaginary.isZero()) {
914             // taking care to preserve the sign of the zero imaginary part
915             return createComplex(FastMath.atan(real), imaginary);
916         } else {
917             // regular formula
918             final FieldComplex<T> n = createComplex(one.add(imaginary), real.negate());
919             final FieldComplex<T> d = createComplex(one.subtract(imaginary),  real);
920             return n.divide(d).log().multiplyPlusI().multiply(0.5);
921         }
922 
923     }
924 
925     /**
926      * Compute the
927      * <a href="http://mathworld.wolfram.com/Cosine.html" TARGET="_top">
928      * cosine</a> of this complex number.
929      * Implements the formula:
930      * <p>
931      *  {@code cos(a + bi) = cos(a)cosh(b) - sin(a)sinh(b)i}
932      * </p><p>
933      * where the (real) functions on the right-hand side are
934      * {@link FastMath#sin}, {@link FastMath#cos},
935      * {@link FastMath#cosh} and {@link FastMath#sinh}.
936      * </p><p>
937      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
938      * input argument is {@code NaN}.
939      * </p><p>
940      * Infinite values in real or imaginary parts of the input may result in
941      * infinite or NaN values returned in parts of the result.</p>
942      * <pre>
943      *  Examples:
944      *  <code>
945      *   cos(1 &plusmn; INFINITY i) = 1 ∓ INFINITY i
946      *   cos(&plusmn;INFINITY + i) = NaN + NaN i
947      *   cos(&plusmn;INFINITY &plusmn; INFINITY i) = NaN + NaN i
948      *  </code>
949      * </pre>
950      *
951      * @return the cosine of this complex number.
952      */
953     @Override
954     public FieldComplex<T> cos() {
955         if (isNaN) {
956             return getNaN(getPartsField());
957         }
958 
959         final FieldSinCos<T>   scr  = FastMath.sinCos(real);
960         final FieldSinhCosh<T> schi = FastMath.sinhCosh(imaginary);
961         return createComplex(scr.cos().multiply(schi.cosh()), scr.sin().negate().multiply(schi.sinh()));
962     }
963 
964     /**
965      * Compute the
966      * <a href="http://mathworld.wolfram.com/HyperbolicCosine.html" TARGET="_top">
967      * hyperbolic cosine</a> of this complex number.
968      * Implements the formula:
969      * <pre>
970      *  <code>
971      *   cosh(a + bi) = cosh(a)cos(b) + sinh(a)sin(b)i
972      *  </code>
973      * </pre>
974      * where the (real) functions on the right-hand side are
975      * {@link FastMath#sin}, {@link FastMath#cos},
976      * {@link FastMath#cosh} and {@link FastMath#sinh}.
977      * <p>
978      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
979      * input argument is {@code NaN}.
980      * </p>
981      * Infinite values in real or imaginary parts of the input may result in
982      * infinite or NaN values returned in parts of the result.
983      * <pre>
984      *  Examples:
985      *  <code>
986      *   cosh(1 &plusmn; INFINITY i) = NaN + NaN i
987      *   cosh(&plusmn;INFINITY + i) = INFINITY &plusmn; INFINITY i
988      *   cosh(&plusmn;INFINITY &plusmn; INFINITY i) = NaN + NaN i
989      *  </code>
990      * </pre>
991      *
992      * @return the hyperbolic cosine of this complex number.
993      */
994     @Override
995     public FieldComplex<T> cosh() {
996         if (isNaN) {
997             return getNaN(getPartsField());
998         }
999 
1000         final FieldSinhCosh<T> schr = FastMath.sinhCosh(real);
1001         final FieldSinCos<T>   sci  = FastMath.sinCos(imaginary);
1002         return createComplex(schr.cosh().multiply(sci.cos()), schr.sinh().multiply(sci.sin()));
1003     }
1004 
1005     /**
1006      * Compute the
1007      * <a href="http://mathworld.wolfram.com/ExponentialFunction.html" TARGET="_top">
1008      * exponential function</a> of this complex number.
1009      * Implements the formula:
1010      * <pre>
1011      *  <code>
1012      *   exp(a + bi) = exp(a)cos(b) + exp(a)sin(b)i
1013      *  </code>
1014      * </pre>
1015      * where the (real) functions on the right-hand side are
1016      * {@link FastMath#exp(CalculusFieldElement)} p}, {@link FastMath#cos(CalculusFieldElement)}, and
1017      * {@link FastMath#sin(CalculusFieldElement)}.
1018      * <p>
1019      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
1020      * input argument is {@code NaN}.
1021      * </p>
1022      * Infinite values in real or imaginary parts of the input may result in
1023      * infinite or NaN values returned in parts of the result.
1024      * <pre>
1025      *  Examples:
1026      *  <code>
1027      *   exp(1 &plusmn; INFINITY i) = NaN + NaN i
1028      *   exp(INFINITY + i) = INFINITY + INFINITY i
1029      *   exp(-INFINITY + i) = 0 + 0i
1030      *   exp(&plusmn;INFINITY &plusmn; INFINITY i) = NaN + NaN i
1031      *  </code>
1032      * </pre>
1033      *
1034      * @return <code><i>e</i><sup>this</sup></code>.
1035      */
1036     @Override
1037     public FieldComplex<T> exp() {
1038         if (isNaN) {
1039             return getNaN(getPartsField());
1040         }
1041 
1042         final T              expReal = FastMath.exp(real);
1043         final FieldSinCos<T> sc      = FastMath.sinCos(imaginary);
1044         return createComplex(expReal.multiply(sc.cos()), expReal.multiply(sc.sin()));
1045     }
1046 
1047     /** {@inheritDoc} */
1048     @Override
1049     public FieldComplex<T> expm1() {
1050         if (isNaN) {
1051             return getNaN(getPartsField());
1052         }
1053 
1054         final T              expm1Real = FastMath.expm1(real);
1055         final FieldSinCos<T> sc        = FastMath.sinCos(imaginary);
1056         return createComplex(expm1Real.multiply(sc.cos()), expm1Real.multiply(sc.sin()));
1057     }
1058 
1059     /**
1060      * Compute the
1061      * <a href="http://mathworld.wolfram.com/NaturalLogarithm.html" TARGET="_top">
1062      * natural logarithm</a> of this complex number.
1063      * Implements the formula:
1064      * <pre>
1065      *  <code>
1066      *   log(a + bi) = ln(|a + bi|) + arg(a + bi)i
1067      *  </code>
1068      * </pre>
1069      * where ln on the right hand side is {@link FastMath#log(CalculusFieldElement)},
1070      * {@code |a + bi|} is the modulus, {@link #abs},  and
1071      * {@code arg(a + bi) = }{@link FastMath#atan2}(b, a).
1072      * <p>
1073      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
1074      * input argument is {@code NaN}.
1075      * </p>
1076      * Infinite (or critical) values in real or imaginary parts of the input may
1077      * result in infinite or NaN values returned in parts of the result.
1078      * <pre>
1079      *  Examples:
1080      *  <code>
1081      *   log(1 &plusmn; INFINITY i) = INFINITY &plusmn; (&pi;/2)i
1082      *   log(INFINITY + i) = INFINITY + 0i
1083      *   log(-INFINITY + i) = INFINITY + &pi;i
1084      *   log(INFINITY &plusmn; INFINITY i) = INFINITY &plusmn; (&pi;/4)i
1085      *   log(-INFINITY &plusmn; INFINITY i) = INFINITY &plusmn; (3&pi;/4)i
1086      *   log(0 + 0i) = -INFINITY + 0i
1087      *  </code>
1088      * </pre>
1089      *
1090      * @return the value <code>ln &nbsp; this</code>, the natural logarithm
1091      * of {@code this}.
1092      */
1093     @Override
1094     public FieldComplex<T> log() {
1095         if (isNaN) {
1096             return getNaN(getPartsField());
1097         }
1098 
1099         return createComplex(FastMath.log(FastMath.hypot(real, imaginary)),
1100                              FastMath.atan2(imaginary, real));
1101     }
1102 
1103     /** {@inheritDoc} */
1104     @Override
1105     public FieldComplex<T> log1p() {
1106         return add(1.0).log();
1107     }
1108 
1109     /** {@inheritDoc} */
1110     @Override
1111     public FieldComplex<T> log10() {
1112         return log().divide(LOG10);
1113     }
1114 
1115     /**
1116      * Returns of value of this complex number raised to the power of {@code x}.
1117      * <p>
1118      * If {@code x} is a real number whose real part has an integer value, returns {@link #pow(int)},
1119      * if both {@code this} and {@code x} are real and {@link FastMath#pow(double, double)}
1120      * with the corresponding real arguments would return a finite number (neither NaN
1121      * nor infinite), then returns the same value converted to {@code Complex},
1122      * with the same special cases.
1123      * In all other cases real cases, implements y<sup>x</sup> = exp(x&middot;log(y)).
1124      * </p>
1125      *
1126      * @param  x exponent to which this {@code Complex} is to be raised.
1127      * @return <code> this<sup>x</sup></code>.
1128      * @throws NullArgumentException if x is {@code null}.
1129      */
1130     @Override
1131     public FieldComplex<T> pow(FieldComplex<T> x)
1132         throws NullArgumentException {
1133 
1134         MathUtils.checkNotNull(x);
1135 
1136         if (x.imaginary.isZero()) {
1137             final int nx = (int) FastMath.rint(x.real.getReal());
1138             if (x.real.getReal() == nx) {
1139                 // integer power
1140                 return pow(nx);
1141             } else if (this.imaginary.isZero()) {
1142                 // check real implementation that handles a bunch of special cases
1143                 final T realPow = FastMath.pow(this.real, x.real);
1144                 if (realPow.isFinite()) {
1145                     return createComplex(realPow, getPartsField().getZero());
1146                 }
1147             }
1148         }
1149 
1150         // generic implementation
1151         return this.log().multiply(x).exp();
1152 
1153     }
1154 
1155 
1156     /**
1157      * Returns of value of this complex number raised to the power of {@code x}.
1158      * <p>
1159      * If {@code x} has an integer value, returns {@link #pow(int)},
1160      * if {@code this} is real and {@link FastMath#pow(double, double)}
1161      * with the corresponding real arguments would return a finite number (neither NaN
1162      * nor infinite), then returns the same value converted to {@code Complex},
1163      * with the same special cases.
1164      * In all other cases real cases, implements y<sup>x</sup> = exp(x&middot;log(y)).
1165      * </p>
1166      *
1167      * @param  x exponent to which this {@code Complex} is to be raised.
1168      * @return <code> this<sup>x</sup></code>.
1169      */
1170     public FieldComplex<T> pow(T x) {
1171 
1172         final int nx = (int) FastMath.rint(x.getReal());
1173         if (x.getReal() == nx) {
1174             // integer power
1175             return pow(nx);
1176         } else if (this.imaginary.isZero()) {
1177             // check real implementation that handles a bunch of special cases
1178             final T realPow = FastMath.pow(this.real, x);
1179             if (realPow.isFinite()) {
1180                 return createComplex(realPow, getPartsField().getZero());
1181             }
1182         }
1183 
1184         // generic implementation
1185         return this.log().multiply(x).exp();
1186 
1187     }
1188 
1189     /**
1190      * Returns of value of this complex number raised to the power of {@code x}.
1191      * <p>
1192      * If {@code x} has an integer value, returns {@link #pow(int)},
1193      * if {@code this} is real and {@link FastMath#pow(double, double)}
1194      * with the corresponding real arguments would return a finite number (neither NaN
1195      * nor infinite), then returns the same value converted to {@code Complex},
1196      * with the same special cases.
1197      * In all other cases real cases, implements y<sup>x</sup> = exp(x&middot;log(y)).
1198      * </p>
1199      *
1200      * @param  x exponent to which this {@code Complex} is to be raised.
1201      * @return <code> this<sup>x</sup></code>.
1202      */
1203     @Override
1204     public FieldComplex<T> pow(double x) {
1205 
1206         final int nx = (int) FastMath.rint(x);
1207         if (x == nx) {
1208             // integer power
1209             return pow(nx);
1210         } else if (this.imaginary.isZero()) {
1211             // check real implementation that handles a bunch of special cases
1212             final T realPow = FastMath.pow(this.real, x);
1213             if (realPow.isFinite()) {
1214                 return createComplex(realPow, getPartsField().getZero());
1215             }
1216         }
1217 
1218         // generic implementation
1219         return this.log().multiply(x).exp();
1220 
1221     }
1222 
1223      /** {@inheritDoc} */
1224     @Override
1225     public FieldComplex<T> pow(final int n) {
1226 
1227         FieldComplex<T> result = getField().getOne();
1228         final boolean invert;
1229         int p = n;
1230         if (p < 0) {
1231             invert = true;
1232             p = -p;
1233         } else {
1234             invert = false;
1235         }
1236 
1237         // Exponentiate by successive squaring
1238         FieldComplex<T> square = this;
1239         while (p > 0) {
1240             if ((p & 0x1) > 0) {
1241                 result = result.multiply(square);
1242             }
1243             square = square.multiply(square);
1244             p = p >> 1;
1245         }
1246 
1247         return invert ? result.reciprocal() : result;
1248 
1249     }
1250 
1251      /**
1252       * Compute the
1253      * <a href="http://mathworld.wolfram.com/Sine.html" TARGET="_top">
1254      * sine</a>
1255      * of this complex number.
1256      * Implements the formula:
1257      * <pre>
1258      *  <code>
1259      *   sin(a + bi) = sin(a)cosh(b) + cos(a)sinh(b)i
1260      *  </code>
1261      * </pre>
1262      * where the (real) functions on the right-hand side are
1263      * {@link FastMath#sin}, {@link FastMath#cos},
1264      * {@link FastMath#cosh} and {@link FastMath#sinh}.
1265      * <p>
1266      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
1267      * input argument is {@code NaN}.
1268      * </p><p>
1269      * Infinite values in real or imaginary parts of the input may result in
1270      * infinite or {@code NaN} values returned in parts of the result.
1271      * <pre>
1272      *  Examples:
1273      *  <code>
1274      *   sin(1 &plusmn; INFINITY i) = 1 &plusmn; INFINITY i
1275      *   sin(&plusmn;INFINITY + i) = NaN + NaN i
1276      *   sin(&plusmn;INFINITY &plusmn; INFINITY i) = NaN + NaN i
1277      *  </code>
1278      * </pre>
1279      *
1280      * @return the sine of this complex number.
1281      */
1282     @Override
1283     public FieldComplex<T> sin() {
1284         if (isNaN) {
1285             return getNaN(getPartsField());
1286         }
1287 
1288         final FieldSinCos<T>   scr  = FastMath.sinCos(real);
1289         final FieldSinhCosh<T> schi = FastMath.sinhCosh(imaginary);
1290         return createComplex(scr.sin().multiply(schi.cosh()), scr.cos().multiply(schi.sinh()));
1291 
1292     }
1293 
1294     /** {@inheritDoc}
1295      */
1296     @Override
1297     public FieldSinCos<FieldComplex<T>> sinCos() {
1298         if (isNaN) {
1299             return new FieldSinCos<>(getNaN(getPartsField()), getNaN(getPartsField()));
1300         }
1301 
1302         final FieldSinCos<T>   scr = FastMath.sinCos(real);
1303         final FieldSinhCosh<T> schi = FastMath.sinhCosh(imaginary);
1304         return new FieldSinCos<>(createComplex(scr.sin().multiply(schi.cosh()), scr.cos().multiply(schi.sinh())),
1305                                  createComplex(scr.cos().multiply(schi.cosh()), scr.sin().negate().multiply(schi.sinh())));
1306     }
1307 
1308     /** {@inheritDoc} */
1309     @Override
1310     public FieldComplex<T> atan2(FieldComplex<T> x) {
1311 
1312         // compute r = sqrt(x^2+y^2)
1313         final FieldComplex<T> r = x.square().add(multiply(this)).sqrt();
1314 
1315         if (x.real.getReal() >= 0) {
1316             // compute atan2(y, x) = 2 atan(y / (r + x))
1317             return divide(r.add(x)).atan().multiply(2);
1318         } else {
1319             // compute atan2(y, x) = +/- pi - 2 atan(y / (r - x))
1320             return divide(r.subtract(x)).atan().multiply(-2).add(x.real.getPi());
1321         }
1322     }
1323 
1324     /** {@inheritDoc}
1325      * <p>
1326      * Branch cuts are on the real axis, below +1.
1327      * </p>
1328      */
1329     @Override
1330     public FieldComplex<T> acosh() {
1331         final FieldComplex<T> sqrtPlus  = add(1).sqrt();
1332         final FieldComplex<T> sqrtMinus = subtract(1).sqrt();
1333         return add(sqrtPlus.multiply(sqrtMinus)).log();
1334     }
1335 
1336     /** {@inheritDoc}
1337      * <p>
1338      * Branch cuts are on the imaginary axis, above +i and below -i.
1339      * </p>
1340      */
1341     @Override
1342     public FieldComplex<T> asinh() {
1343         return add(multiply(this).add(1.0).sqrt()).log();
1344     }
1345 
1346     /** {@inheritDoc}
1347      * <p>
1348      * Branch cuts are on the real axis, above +1 and below -1.
1349      * </p>
1350      */
1351     @Override
1352     public FieldComplex<T> atanh() {
1353         final FieldComplex<T> logPlus  = add(1).log();
1354         final FieldComplex<T> logMinus = createComplex(getPartsField().getOne().subtract(real), imaginary.negate()).log();
1355         return logPlus.subtract(logMinus).multiply(0.5);
1356     }
1357 
1358     /**
1359      * Compute the
1360      * <a href="http://mathworld.wolfram.com/HyperbolicSine.html" TARGET="_top">
1361      * hyperbolic sine</a> of this complex number.
1362      * Implements the formula:
1363      * <pre>
1364      *  <code>
1365      *   sinh(a + bi) = sinh(a)cos(b)) + cosh(a)sin(b)i
1366      *  </code>
1367      * </pre>
1368      * where the (real) functions on the right-hand side are
1369      * {@link FastMath#sin}, {@link FastMath#cos},
1370      * {@link FastMath#cosh} and {@link FastMath#sinh}.
1371      * <p>
1372      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
1373      * input argument is {@code NaN}.
1374      * </p><p>
1375      * Infinite values in real or imaginary parts of the input may result in
1376      * infinite or NaN values returned in parts of the result.
1377      * <pre>
1378      *  Examples:
1379      *  <code>
1380      *   sinh(1 &plusmn; INFINITY i) = NaN + NaN i
1381      *   sinh(&plusmn;INFINITY + i) = &plusmn; INFINITY + INFINITY i
1382      *   sinh(&plusmn;INFINITY &plusmn; INFINITY i) = NaN + NaN i
1383      *  </code>
1384      * </pre>
1385      *
1386      * @return the hyperbolic sine of {@code this}.
1387      */
1388     @Override
1389     public FieldComplex<T> sinh() {
1390         if (isNaN) {
1391             return getNaN(getPartsField());
1392         }
1393 
1394         final FieldSinhCosh<T> schr = FastMath.sinhCosh(real);
1395         final FieldSinCos<T>   sci  = FastMath.sinCos(imaginary);
1396         return createComplex(schr.sinh().multiply(sci.cos()), schr.cosh().multiply(sci.sin()));
1397     }
1398 
1399     /** {@inheritDoc}
1400      */
1401     @Override
1402     public FieldSinhCosh<FieldComplex<T>> sinhCosh() {
1403         if (isNaN) {
1404             return new FieldSinhCosh<>(getNaN(getPartsField()), getNaN(getPartsField()));
1405         }
1406 
1407         final FieldSinhCosh<T> schr = FastMath.sinhCosh(real);
1408         final FieldSinCos<T>   sci  = FastMath.sinCos(imaginary);
1409         return new FieldSinhCosh<>(createComplex(schr.sinh().multiply(sci.cos()), schr.cosh().multiply(sci.sin())),
1410                                    createComplex(schr.cosh().multiply(sci.cos()), schr.sinh().multiply(sci.sin())));
1411     }
1412 
1413     /**
1414      * Compute the
1415      * <a href="http://mathworld.wolfram.com/SquareRoot.html" TARGET="_top">
1416      * square root</a> of this complex number.
1417      * Implements the following algorithm to compute {@code sqrt(a + bi)}:
1418      * <ol><li>Let {@code t = sqrt((|a| + |a + bi|) / 2)}</li>
1419      * <li><pre>if {@code  a ≥ 0} return {@code t + (b/2t)i}
1420      *  else return {@code |b|/2t + sign(b)t i }</pre></li>
1421      * </ol>
1422      * where <ul>
1423      * <li>{@code |a| = }{@link FastMath#abs(CalculusFieldElement) abs(a)}</li>
1424      * <li>{@code |a + bi| = }{@link FastMath#hypot(CalculusFieldElement, CalculusFieldElement) hypot(a, b)}</li>
1425      * <li>{@code sign(b) = }{@link FastMath#copySign(CalculusFieldElement, CalculusFieldElement) copySign(1, b)}
1426      * </ul>
1427      * The real part is therefore always nonnegative.
1428      * <p>
1429      * Returns {@link #getNaN(Field) NaN} if either real or imaginary part of the
1430      * input argument is {@code NaN}.
1431      * </p>
1432      * <p>
1433      * Infinite values in real or imaginary parts of the input may result in
1434      * infinite or NaN values returned in parts of the result.
1435      * </p>
1436      * <pre>
1437      *  Examples:
1438      *  <code>
1439      *   sqrt(1 ± ∞ i) = ∞ + NaN i
1440      *   sqrt(∞ + i) = ∞ + 0i
1441      *   sqrt(-∞ + i) = 0 + ∞ i
1442      *   sqrt(∞ ± ∞ i) = ∞ + NaN i
1443      *   sqrt(-∞ ± ∞ i) = NaN ± ∞ i
1444      *  </code>
1445      * </pre>
1446      *
1447      * @return the square root of {@code this} with nonnegative real part.
1448      */
1449     @Override
1450     public FieldComplex<T> sqrt() {
1451         if (isNaN) {
1452             return getNaN(getPartsField());
1453         }
1454 
1455         if (isZero()) {
1456             return getZero(getPartsField());
1457         }
1458 
1459         T t = FastMath.sqrt((FastMath.abs(real).add(FastMath.hypot(real, imaginary))).multiply(0.5));
1460         if (real.getReal() >= 0.0) {
1461             return createComplex(t, imaginary.divide(t.multiply(2)));
1462         } else {
1463             return createComplex(FastMath.abs(imaginary).divide(t.multiply(2)),
1464                                  FastMath.copySign(t, imaginary));
1465         }
1466     }
1467 
1468     /**
1469      * Compute the
1470      * <a href="http://mathworld.wolfram.com/SquareRoot.html" TARGET="_top">
1471      * square root</a> of <code>1 - this<sup>2</sup></code> for this complex
1472      * number.
1473      * Computes the result directly as
1474      * {@code sqrt(ONE.subtract(z.square()))}.
1475      * <p>
1476      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
1477      * input argument is {@code NaN}.
1478      * </p>
1479      * Infinite values in real or imaginary parts of the input may result in
1480      * infinite or NaN values returned in parts of the result.
1481      *
1482      * @return the square root of <code>1 - this<sup>2</sup></code>.
1483      */
1484     public FieldComplex<T> sqrt1z() {
1485         final FieldComplex<T> t2 = this.square();
1486         return createComplex(getPartsField().getOne().subtract(t2.real), t2.imaginary.negate()).sqrt();
1487     }
1488 
1489     /** {@inheritDoc}
1490      * <p>
1491      * This implementation compute the principal cube root by using a branch cut along real negative axis.
1492      * </p>
1493      */
1494     @Override
1495     public FieldComplex<T> cbrt() {
1496         final T              magnitude = FastMath.cbrt(abs().getRealPart());
1497         final FieldSinCos<T> sc        = FastMath.sinCos(getArgument().divide(3));
1498         return createComplex(magnitude.multiply(sc.cos()), magnitude.multiply(sc.sin()));
1499     }
1500 
1501     /** {@inheritDoc}
1502      * <p>
1503      * This implementation compute the principal n<sup>th</sup> root by using a branch cut along real negative axis.
1504      * </p>
1505      */
1506     @Override
1507     public FieldComplex<T> rootN(int n) {
1508         final T              magnitude = FastMath.pow(abs().getRealPart(), 1.0 / n);
1509         final FieldSinCos<T> sc        = FastMath.sinCos(getArgument().divide(n));
1510         return createComplex(magnitude.multiply(sc.cos()), magnitude.multiply(sc.sin()));
1511     }
1512 
1513     /**
1514      * Compute the
1515      * <a href="http://mathworld.wolfram.com/Tangent.html" TARGET="_top">
1516      * tangent</a> of this complex number.
1517      * Implements the formula:
1518      * <pre>
1519      *  <code>
1520      *   tan(a + bi) = sin(2a)/(cos(2a)+cosh(2b)) + [sinh(2b)/(cos(2a)+cosh(2b))]i
1521      *  </code>
1522      * </pre>
1523      * where the (real) functions on the right-hand side are
1524      * {@link FastMath#sin}, {@link FastMath#cos}, {@link FastMath#cosh} and
1525      * {@link FastMath#sinh}.
1526      * <p>
1527      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
1528      * input argument is {@code NaN}.
1529      * </p>
1530      * Infinite (or critical) values in real or imaginary parts of the input may
1531      * result in infinite or NaN values returned in parts of the result.
1532      * <pre>
1533      *  Examples:
1534      *  <code>
1535      *   tan(a &plusmn; INFINITY i) = 0 &plusmn; i
1536      *   tan(&plusmn;INFINITY + bi) = NaN + NaN i
1537      *   tan(&plusmn;INFINITY &plusmn; INFINITY i) = NaN + NaN i
1538      *   tan(&plusmn;&pi;/2 + 0 i) = &plusmn;INFINITY + NaN i
1539      *  </code>
1540      * </pre>
1541      *
1542      * @return the tangent of {@code this}.
1543      */
1544     @Override
1545     public FieldComplex<T> tan() {
1546         if (isNaN || real.isInfinite()) {
1547             return getNaN(getPartsField());
1548         }
1549         if (imaginary.getReal() > 20.0) {
1550             return getI(getPartsField());
1551         }
1552         if (imaginary.getReal() < -20.0) {
1553             return getMinusI(getPartsField());
1554         }
1555 
1556         final FieldSinCos<T> sc2r = FastMath.sinCos(real.multiply(2));
1557         T imaginary2 = imaginary.multiply(2);
1558         T d = sc2r.cos().add(FastMath.cosh(imaginary2));
1559 
1560         return createComplex(sc2r.sin().divide(d), FastMath.sinh(imaginary2).divide(d));
1561 
1562     }
1563 
1564     /**
1565      * Compute the
1566      * <a href="http://mathworld.wolfram.com/HyperbolicTangent.html" TARGET="_top">
1567      * hyperbolic tangent</a> of this complex number.
1568      * Implements the formula:
1569      * <pre>
1570      *  <code>
1571      *   tan(a + bi) = sinh(2a)/(cosh(2a)+cos(2b)) + [sin(2b)/(cosh(2a)+cos(2b))]i
1572      *  </code>
1573      * </pre>
1574      * where the (real) functions on the right-hand side are
1575      * {@link FastMath#sin}, {@link FastMath#cos}, {@link FastMath#cosh} and
1576      * {@link FastMath#sinh}.
1577      * <p>
1578      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
1579      * input argument is {@code NaN}.
1580      * </p>
1581      * Infinite values in real or imaginary parts of the input may result in
1582      * infinite or NaN values returned in parts of the result.
1583      * <pre>
1584      *  Examples:
1585      *  <code>
1586      *   tanh(a &plusmn; INFINITY i) = NaN + NaN i
1587      *   tanh(&plusmn;INFINITY + bi) = &plusmn;1 + 0 i
1588      *   tanh(&plusmn;INFINITY &plusmn; INFINITY i) = NaN + NaN i
1589      *   tanh(0 + (&pi;/2)i) = NaN + INFINITY i
1590      *  </code>
1591      * </pre>
1592      *
1593      * @return the hyperbolic tangent of {@code this}.
1594      */
1595     @Override
1596     public FieldComplex<T> tanh() {
1597         if (isNaN || imaginary.isInfinite()) {
1598             return getNaN(getPartsField());
1599         }
1600         if (real.getReal() > 20.0) {
1601             return getOne(getPartsField());
1602         }
1603         if (real.getReal() < -20.0) {
1604             return getMinusOne(getPartsField());
1605         }
1606         T real2 = real.multiply(2);
1607         final FieldSinCos<T> sc2i = FastMath.sinCos(imaginary.multiply(2));
1608         T d = FastMath.cosh(real2).add(sc2i.cos());
1609 
1610         return createComplex(FastMath.sinh(real2).divide(d), sc2i.sin().divide(d));
1611     }
1612 
1613 
1614 
1615     /**
1616      * Compute the argument of this complex number.
1617      * The argument is the angle phi between the positive real axis and
1618      * the point representing this number in the complex plane.
1619      * The value returned is between -PI (not inclusive)
1620      * and PI (inclusive), with negative values returned for numbers with
1621      * negative imaginary parts.
1622      * <p>
1623      * If either real or imaginary part (or both) is NaN, NaN is returned.
1624      * Infinite parts are handled as {@code Math.atan2} handles them,
1625      * essentially treating finite parts as zero in the presence of an
1626      * infinite coordinate and returning a multiple of pi/4 depending on
1627      * the signs of the infinite parts.
1628      * See the javadoc for {@code Math.atan2} for full details.
1629      *
1630      * @return the argument of {@code this}.
1631      */
1632     public T getArgument() {
1633         return FastMath.atan2(getImaginaryPart(), getRealPart());
1634     }
1635 
1636     /**
1637      * Computes the n-th roots of this complex number.
1638      * The nth roots are defined by the formula:
1639      * <pre>
1640      *  <code>
1641      *   z<sub>k</sub> = abs<sup>1/n</sup> (cos(phi + 2&pi;k/n) + i (sin(phi + 2&pi;k/n))
1642      *  </code>
1643      * </pre>
1644      * for <i>{@code k=0, 1, ..., n-1}</i>, where {@code abs} and {@code phi}
1645      * are respectively the {@link #abs() modulus} and
1646      * {@link #getArgument() argument} of this complex number.
1647      * <p>
1648      * If one or both parts of this complex number is NaN, a list with just
1649      * one element, {@link #getNaN(Field)} is returned.
1650      * if neither part is NaN, but at least one part is infinite, the result
1651      * is a one-element list containing {@link #getInf(Field)}.
1652      *
1653      * @param n Degree of root.
1654      * @return a List of all {@code n}-th roots of {@code this}.
1655      * @throws MathIllegalArgumentException if {@code n <= 0}.
1656      */
1657     public List<FieldComplex<T>> nthRoot(int n) throws MathIllegalArgumentException {
1658 
1659         if (n <= 0) {
1660             throw new MathIllegalArgumentException(LocalizedCoreFormats.CANNOT_COMPUTE_NTH_ROOT_FOR_NEGATIVE_N,
1661                                                    n);
1662         }
1663 
1664         final List<FieldComplex<T>> result = new ArrayList<>();
1665 
1666         if (isNaN) {
1667             result.add(getNaN(getPartsField()));
1668             return result;
1669         }
1670         if (isInfinite()) {
1671             result.add(getInf(getPartsField()));
1672             return result;
1673         }
1674 
1675         // nth root of abs -- faster / more accurate to use a solver here?
1676         final T nthRootOfAbs = FastMath.pow(FastMath.hypot(real, imaginary), 1.0 / n);
1677 
1678         // Compute nth roots of complex number with k = 0, 1, ... n-1
1679         final T nthPhi = getArgument().divide(n);
1680         final double slice = 2 * FastMath.PI / n;
1681         T innerPart = nthPhi;
1682         for (int k = 0; k < n ; k++) {
1683             // inner part
1684             final FieldSinCos<T> scInner = FastMath.sinCos(innerPart);
1685             final T realPart = nthRootOfAbs.multiply(scInner.cos());
1686             final T imaginaryPart = nthRootOfAbs.multiply(scInner.sin());
1687             result.add(createComplex(realPart, imaginaryPart));
1688             innerPart = innerPart.add(slice);
1689         }
1690 
1691         return result;
1692     }
1693 
1694     /**
1695      * Create a complex number given the real and imaginary parts.
1696      *
1697      * @param realPart Real part.
1698      * @param imaginaryPart Imaginary part.
1699      * @return a new complex number instance.
1700      *
1701      * @see #valueOf(CalculusFieldElement, CalculusFieldElement)
1702      */
1703     protected FieldComplex<T> createComplex(final T realPart, final T imaginaryPart) {
1704         return new FieldComplex<>(realPart, imaginaryPart);
1705     }
1706 
1707     /**
1708      * Create a complex number given the real and imaginary parts.
1709      *
1710      * @param realPart Real part.
1711      * @param imaginaryPart Imaginary part.
1712      * @param <T> the type of the field elements
1713      * @return a Complex instance.
1714      */
1715     public static <T extends CalculusFieldElement<T>> FieldComplex<T>
1716         valueOf(T realPart, T imaginaryPart) {
1717         if (realPart.isNaN() || imaginaryPart.isNaN()) {
1718             return getNaN(realPart.getField());
1719         }
1720         return new FieldComplex<>(realPart, imaginaryPart);
1721     }
1722 
1723     /**
1724      * Create a complex number given only the real part.
1725      *
1726      * @param realPart Real part.
1727      * @param <T> the type of the field elements
1728      * @return a Complex instance.
1729      */
1730     public static <T extends CalculusFieldElement<T>> FieldComplex<T>
1731         valueOf(T realPart) {
1732         if (realPart.isNaN()) {
1733             return getNaN(realPart.getField());
1734         }
1735         return new FieldComplex<>(realPart);
1736     }
1737 
1738     /** {@inheritDoc} */
1739     @Override
1740     public FieldComplex<T> newInstance(double realPart) {
1741         return valueOf(getPartsField().getZero().newInstance(realPart));
1742     }
1743 
1744     /** {@inheritDoc} */
1745     @Override
1746     public FieldComplexField<T> getField() {
1747         return FieldComplexField.getField(getPartsField());
1748     }
1749 
1750     /** Get the {@link Field} the real and imaginary parts belong to.
1751      * @return {@link Field} the real and imaginary parts belong to
1752      */
1753     public Field<T> getPartsField() {
1754         return real.getField();
1755     }
1756 
1757     /** {@inheritDoc} */
1758     @Override
1759     public String toString() {
1760         return "(" + real + ", " + imaginary + ")";
1761     }
1762 
1763     /** {@inheritDoc} */
1764     @Override
1765     public FieldComplex<T> scalb(int n) {
1766         return createComplex(FastMath.scalb(real, n), FastMath.scalb(imaginary, n));
1767     }
1768 
1769     /** {@inheritDoc} */
1770     @Override
1771     public FieldComplex<T> ulp() {
1772         return createComplex(FastMath.ulp(real), FastMath.ulp(imaginary));
1773     }
1774 
1775     /** {@inheritDoc} */
1776     @Override
1777     public FieldComplex<T> hypot(FieldComplex<T> y) {
1778         if (isInfinite() || y.isInfinite()) {
1779             return getInf(getPartsField());
1780         } else if (isNaN() || y.isNaN()) {
1781             return getNaN(getPartsField());
1782         } else {
1783             return square().add(y.square()).sqrt();
1784         }
1785     }
1786 
1787     /** {@inheritDoc} */
1788     @Override
1789     public FieldComplex<T> linearCombination(final FieldComplex<T>[] a, final FieldComplex<T>[] b)
1790         throws MathIllegalArgumentException {
1791         final int n = 2 * a.length;
1792         final T[] realA      = MathArrays.buildArray(getPartsField(), n);
1793         final T[] realB      = MathArrays.buildArray(getPartsField(), n);
1794         final T[] imaginaryA = MathArrays.buildArray(getPartsField(), n);
1795         final T[] imaginaryB = MathArrays.buildArray(getPartsField(), n);
1796         for (int i = 0; i < a.length; ++i)  {
1797             final FieldComplex<T> ai = a[i];
1798             final FieldComplex<T> bi = b[i];
1799             realA[2 * i    ]      = ai.real;
1800             realA[2 * i + 1]      = ai.imaginary.negate();
1801             realB[2 * i    ]      = bi.real;
1802             realB[2 * i + 1]      = bi.imaginary;
1803             imaginaryA[2 * i    ] = ai.real;
1804             imaginaryA[2 * i + 1] = ai.imaginary;
1805             imaginaryB[2 * i    ] = bi.imaginary;
1806             imaginaryB[2 * i + 1] = bi.real;
1807         }
1808         return createComplex(real.linearCombination(realA,  realB),
1809                              real.linearCombination(imaginaryA, imaginaryB));
1810     }
1811 
1812     /** {@inheritDoc} */
1813     @Override
1814     public FieldComplex<T> linearCombination(final double[] a, final FieldComplex<T>[] b)
1815         throws MathIllegalArgumentException {
1816         final int n = a.length;
1817         final T[] realB      = MathArrays.buildArray(getPartsField(), n);
1818         final T[] imaginaryB = MathArrays.buildArray(getPartsField(), n);
1819         for (int i = 0; i < a.length; ++i)  {
1820             final FieldComplex<T> bi = b[i];
1821             realB[i]      = bi.real;
1822             imaginaryB[i] = bi.imaginary;
1823         }
1824         return createComplex(real.linearCombination(a,  realB),
1825                              real.linearCombination(a, imaginaryB));
1826     }
1827 
1828     /** {@inheritDoc} */
1829     @Override
1830     public FieldComplex<T> linearCombination(final FieldComplex<T> a1, final FieldComplex<T> b1, final FieldComplex<T> a2, final FieldComplex<T> b2) {
1831         return createComplex(real.linearCombination(a1.real, b1.real,
1832                                                     a1.imaginary.negate(), b1.imaginary,
1833                                                     a2.real, b2.real,
1834                                                     a2.imaginary.negate(), b2.imaginary),
1835                              real.linearCombination(a1.real, b1.imaginary,
1836                                                     a1.imaginary, b1.real,
1837                                                     a2.real, b2.imaginary,
1838                                                     a2.imaginary, b2.real));
1839     }
1840 
1841     /** {@inheritDoc} */
1842     @Override
1843     public FieldComplex<T> linearCombination(final double a1, final FieldComplex<T> b1, final double a2, final FieldComplex<T> b2) {
1844         return createComplex(real.linearCombination(a1, b1.real,
1845                                                     a2, b2.real),
1846                              real.linearCombination(a1, b1.imaginary,
1847                                                     a2, b2.imaginary));
1848     }
1849 
1850     /** {@inheritDoc} */
1851     @Override
1852     public FieldComplex<T> linearCombination(final FieldComplex<T> a1, final FieldComplex<T> b1,
1853                                                 final FieldComplex<T> a2, final FieldComplex<T> b2,
1854                                                 final FieldComplex<T> a3, final FieldComplex<T> b3) {
1855         FieldComplex<T>[] a = MathArrays.buildArray(getField(), 3);
1856         a[0] = a1;
1857         a[1] = a2;
1858         a[2] = a3;
1859         FieldComplex<T>[] b = MathArrays.buildArray(getField(), 3);
1860         b[0] = b1;
1861         b[1] = b2;
1862         b[2] = b3;
1863         return linearCombination(a, b);
1864     }
1865 
1866     /** {@inheritDoc} */
1867     @Override
1868     public FieldComplex<T> linearCombination(final double a1, final FieldComplex<T> b1,
1869                                                 final double a2, final FieldComplex<T> b2,
1870                                                 final double a3, final FieldComplex<T> b3) {
1871         FieldComplex<T>[] b = MathArrays.buildArray(getField(), 3);
1872         b[0] = b1;
1873         b[1] = b2;
1874         b[2] = b3;
1875         return linearCombination(new double[]  { a1, a2, a3 }, b);
1876     }
1877 
1878     /** {@inheritDoc} */
1879     @Override
1880     public FieldComplex<T> linearCombination(final FieldComplex<T> a1, final FieldComplex<T> b1,
1881                                                 final FieldComplex<T> a2, final FieldComplex<T> b2,
1882                                                 final FieldComplex<T> a3, final FieldComplex<T> b3,
1883                                                 final FieldComplex<T> a4, final FieldComplex<T> b4) {
1884         FieldComplex<T>[] a = MathArrays.buildArray(getField(), 4);
1885         a[0] = a1;
1886         a[1] = a2;
1887         a[2] = a3;
1888         a[3] = a4;
1889         FieldComplex<T>[] b = MathArrays.buildArray(getField(), 4);
1890         b[0] = b1;
1891         b[1] = b2;
1892         b[2] = b3;
1893         b[3] = b4;
1894         return linearCombination(a, b);
1895     }
1896 
1897     /** {@inheritDoc} */
1898     @Override
1899     public FieldComplex<T> linearCombination(final double a1, final FieldComplex<T> b1,
1900                                                 final double a2, final FieldComplex<T> b2,
1901                                                 final double a3, final FieldComplex<T> b3,
1902                                                 final double a4, final FieldComplex<T> b4) {
1903         FieldComplex<T>[] b = MathArrays.buildArray(getField(), 4);
1904         b[0] = b1;
1905         b[1] = b2;
1906         b[2] = b3;
1907         b[3] = b4;
1908         return linearCombination(new double[]  { a1, a2, a3, a4 }, b);
1909     }
1910 
1911     /** {@inheritDoc} */
1912     @Override
1913     public FieldComplex<T> ceil() {
1914         return createComplex(FastMath.ceil(getRealPart()), FastMath.ceil(getImaginaryPart()));
1915     }
1916 
1917     /** {@inheritDoc} */
1918     @Override
1919     public FieldComplex<T> floor() {
1920         return createComplex(FastMath.floor(getRealPart()), FastMath.floor(getImaginaryPart()));
1921     }
1922 
1923     /** {@inheritDoc} */
1924     @Override
1925     public FieldComplex<T> rint() {
1926         return createComplex(FastMath.rint(getRealPart()), FastMath.rint(getImaginaryPart()));
1927     }
1928 
1929     /** {@inheritDoc}
1930      * <p>
1931      * for complex numbers, the integer n corresponding to {@code this.subtract(remainder(a)).divide(a)}
1932      * is a <a href="https://en.wikipedia.org/wiki/Gaussian_integer">Wikipedia - Gaussian integer</a>.
1933      * </p>
1934      */
1935     @Override
1936     public FieldComplex<T> remainder(final double a) {
1937         return createComplex(FastMath.IEEEremainder(getRealPart(), a), FastMath.IEEEremainder(getImaginaryPart(), a));
1938     }
1939 
1940     /** {@inheritDoc}
1941      * <p>
1942      * for complex numbers, the integer n corresponding to {@code this.subtract(remainder(a)).divide(a)}
1943      * is a <a href="https://en.wikipedia.org/wiki/Gaussian_integer">Wikipedia - Gaussian integer</a>.
1944      * </p>
1945      */
1946     @Override
1947     public FieldComplex<T> remainder(final FieldComplex<T> a) {
1948         final FieldComplex<T> complexQuotient = divide(a);
1949         final T  qRInt           = FastMath.rint(complexQuotient.real);
1950         final T  qIInt           = FastMath.rint(complexQuotient.imaginary);
1951         return createComplex(real.subtract(qRInt.multiply(a.real)).add(qIInt.multiply(a.imaginary)),
1952                              imaginary.subtract(qRInt.multiply(a.imaginary)).subtract(qIInt.multiply(a.real)));
1953     }
1954 
1955     /** {@inheritDoc} */
1956     @Override
1957     public FieldComplex<T> sign() {
1958         if (isNaN() || isZero()) {
1959             return this;
1960         } else {
1961             return this.divide(FastMath.hypot(real, imaginary));
1962         }
1963     }
1964 
1965     /** {@inheritDoc}
1966      * <p>
1967      * The signs of real and imaginary parts are copied independently.
1968      * </p>
1969      */
1970     @Override
1971     public FieldComplex<T> copySign(final FieldComplex<T> z) {
1972         return createComplex(FastMath.copySign(getRealPart(), z.getRealPart()),
1973                              FastMath.copySign(getImaginaryPart(), z.getImaginaryPart()));
1974     }
1975 
1976     /** {@inheritDoc} */
1977     @Override
1978     public FieldComplex<T> copySign(double r) {
1979         return createComplex(FastMath.copySign(getRealPart(), r), FastMath.copySign(getImaginaryPart(), r));
1980     }
1981 
1982     /** {@inheritDoc} */
1983     @Override
1984     public FieldComplex<T> toDegrees() {
1985         return createComplex(FastMath.toDegrees(getRealPart()), FastMath.toDegrees(getImaginaryPart()));
1986     }
1987 
1988     /** {@inheritDoc} */
1989     @Override
1990     public FieldComplex<T> toRadians() {
1991         return createComplex(FastMath.toRadians(getRealPart()), FastMath.toRadians(getImaginaryPart()));
1992     }
1993 
1994     /** {@inheritDoc} */
1995     @Override
1996     public FieldComplex<T> getPi() {
1997         return getPi(getPartsField());
1998     }
1999 
2000 }