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     /**
584      * Access the real part.
585      *
586      * @return the real part.
587      */
588     public T getRealPart() {
589         return real;
590     }
591 
592     /**
593      * Checks whether either or both parts of this complex number is
594      * {@code NaN}.
595      *
596      * @return true if either or both parts of this complex number is
597      * {@code NaN}; false otherwise.
598      */
599     @Override
600     public boolean isNaN() {
601         return isNaN;
602     }
603 
604     /** Check whether the instance is real (i.e. imaginary part is zero).
605      * @return true if imaginary part is zero
606       */
607     public boolean isReal() {
608         return imaginary.isZero();
609     }
610 
611     /** Check whether the instance is an integer (i.e. imaginary part is zero and real part has no fractional part).
612      * @return true if imaginary part is zero and real part has no fractional part
613      */
614     public boolean isMathematicalInteger() {
615         return isReal() && Precision.isMathematicalInteger(real.getReal());
616     }
617 
618     /**
619      * Checks whether either the real or imaginary part of this complex number
620      * takes an infinite value (either {@code Double.POSITIVE_INFINITY} or
621      * {@code Double.NEGATIVE_INFINITY}) and neither part
622      * is {@code NaN}.
623      *
624      * @return true if one or both parts of this complex number are infinite
625      * and neither part is {@code NaN}.
626      */
627     @Override
628     public boolean isInfinite() {
629         return isInfinite;
630     }
631 
632     /**
633      * Returns a {@code Complex} whose value is {@code this * factor}.
634      * Implements preliminary checks for {@code NaN} and infinity followed by
635      * the definitional formula:
636      * <p>
637      *   {@code (a + bi)(c + di) = (ac - bd) + (ad + bc)i}
638      * </p>
639      * Returns {@link #getNaN(Field)} if either {@code this} or {@code factor} has one or
640      * more {@code NaN} parts.
641      * <p>
642      * Returns {@link #getInf(Field)} if neither {@code this} nor {@code factor} has one
643      * or more {@code NaN} parts and if either {@code this} or {@code factor}
644      * has one or more infinite parts (same result is returned regardless of
645      * the sign of the components).
646      * </p><p>
647      * Returns finite values in components of the result per the definitional
648      * formula in all remaining cases.</p>
649      *
650      * @param  factor value to be multiplied by this {@code Complex}.
651      * @return {@code this * factor}.
652      * @throws NullArgumentException if {@code factor} is {@code null}.
653      */
654     @Override
655     public FieldComplex<T> multiply(FieldComplex<T> factor)
656         throws NullArgumentException {
657         MathUtils.checkNotNull(factor);
658         if (isNaN || factor.isNaN) {
659             return getNaN(getPartsField());
660         }
661         if (real.isInfinite() ||
662             imaginary.isInfinite() ||
663             factor.real.isInfinite() ||
664             factor.imaginary.isInfinite()) {
665             // we don't use isInfinite() to avoid testing for NaN again
666             return getInf(getPartsField());
667         }
668         return createComplex(real.linearCombination(real, factor.real, imaginary.negate(), factor.imaginary),
669                              real.linearCombination(real, factor.imaginary, imaginary, factor.real));
670     }
671 
672     /**
673      * Returns a {@code Complex} whose value is {@code this * factor}, with {@code factor}
674      * interpreted as a integer number.
675      *
676      * @param  factor value to be multiplied by this {@code Complex}.
677      * @return {@code this * factor}.
678      * @see #multiply(FieldComplex)
679      */
680     @Override
681     public FieldComplex<T> multiply(final int factor) {
682         if (isNaN) {
683             return getNaN(getPartsField());
684         }
685         if (real.isInfinite() || imaginary.isInfinite()) {
686             return getInf(getPartsField());
687         }
688         return createComplex(real.multiply(factor), imaginary.multiply(factor));
689     }
690 
691     /**
692      * Returns a {@code Complex} whose value is {@code this * factor}, with {@code factor}
693      * interpreted as a real number.
694      *
695      * @param  factor value to be multiplied by this {@code Complex}.
696      * @return {@code this * factor}.
697      * @see #multiply(FieldComplex)
698      */
699     @Override
700     public FieldComplex<T> multiply(double factor) {
701         if (isNaN || Double.isNaN(factor)) {
702             return getNaN(getPartsField());
703         }
704         if (real.isInfinite() ||
705             imaginary.isInfinite() ||
706             Double.isInfinite(factor)) {
707             // we don't use isInfinite() to avoid testing for NaN again
708             return getInf(getPartsField());
709         }
710         return createComplex(real.multiply(factor), imaginary.multiply(factor));
711     }
712 
713     /**
714      * Returns a {@code Complex} whose value is {@code this * factor}, with {@code factor}
715      * interpreted as a real number.
716      *
717      * @param  factor value to be multiplied by this {@code Complex}.
718      * @return {@code this * factor}.
719      * @see #multiply(FieldComplex)
720      */
721     public FieldComplex<T> multiply(T factor) {
722         if (isNaN || factor.isNaN()) {
723             return getNaN(getPartsField());
724         }
725         if (real.isInfinite() ||
726             imaginary.isInfinite() ||
727             factor.isInfinite()) {
728             // we don't use isInfinite() to avoid testing for NaN again
729             return getInf(getPartsField());
730         }
731         return createComplex(real.multiply(factor), imaginary.multiply(factor));
732     }
733 
734     /** Compute this * i.
735      * @return this * i
736      * @since 2.0
737      */
738     public FieldComplex<T> multiplyPlusI() {
739         return createComplex(imaginary.negate(), real);
740     }
741 
742     /** Compute this *- -i.
743      * @return this * i
744      * @since 2.0
745      */
746     public FieldComplex<T> multiplyMinusI() {
747         return createComplex(imaginary, real.negate());
748     }
749 
750     @Override
751     public FieldComplex<T> square() {
752         return multiply(this);
753     }
754 
755     /**
756      * Returns a {@code Complex} whose value is {@code (-this)}.
757      * Returns {@code NaN} if either real or imaginary
758      * part of this Complex number is {@code Double.NaN}.
759      *
760      * @return {@code -this}.
761      */
762     @Override
763     public FieldComplex<T> negate() {
764         if (isNaN) {
765             return getNaN(getPartsField());
766         }
767 
768         return createComplex(real.negate(), imaginary.negate());
769     }
770 
771     /**
772      * Returns a {@code Complex} whose value is
773      * {@code (this - subtrahend)}.
774      * Uses the definitional formula
775      * <p>
776      *  {@code (a + bi) - (c + di) = (a-c) + (b-d)i}
777      * </p>
778      * If either {@code this} or {@code subtrahend} has a {@code NaN]} value in either part,
779      * {@link #getNaN(Field)} is returned; otherwise infinite and {@code NaN} values are
780      * returned in the parts of the result according to the rules for
781      * {@link java.lang.Double} arithmetic.
782      *
783      * @param  subtrahend value to be subtracted from this {@code Complex}.
784      * @return {@code this - subtrahend}.
785      * @throws NullArgumentException if {@code subtrahend} is {@code null}.
786      */
787     @Override
788     public FieldComplex<T> subtract(FieldComplex<T> subtrahend)
789         throws NullArgumentException {
790         MathUtils.checkNotNull(subtrahend);
791         if (isNaN || subtrahend.isNaN) {
792             return getNaN(getPartsField());
793         }
794 
795         return createComplex(real.subtract(subtrahend.getRealPart()),
796                              imaginary.subtract(subtrahend.getImaginaryPart()));
797     }
798 
799     /**
800      * Returns a {@code Complex} whose value is
801      * {@code (this - subtrahend)}.
802      *
803      * @param  subtrahend value to be subtracted from this {@code Complex}.
804      * @return {@code this - subtrahend}.
805      * @see #subtract(FieldComplex)
806      */
807     @Override
808     public FieldComplex<T> subtract(double subtrahend) {
809         if (isNaN || Double.isNaN(subtrahend)) {
810             return getNaN(getPartsField());
811         }
812         return createComplex(real.subtract(subtrahend), imaginary);
813     }
814 
815     /**
816      * Returns a {@code Complex} whose value is
817      * {@code (this - subtrahend)}.
818      *
819      * @param  subtrahend value to be subtracted from this {@code Complex}.
820      * @return {@code this - subtrahend}.
821      * @see #subtract(FieldComplex)
822      */
823     public FieldComplex<T> subtract(T subtrahend) {
824         if (isNaN || subtrahend.isNaN()) {
825             return getNaN(getPartsField());
826         }
827         return createComplex(real.subtract(subtrahend), imaginary);
828     }
829 
830     /**
831      * Compute the
832      * <a href="http://mathworld.wolfram.com/InverseCosine.html" TARGET="_top">
833      * inverse cosine</a> of this complex number.
834      * Implements the formula:
835      * <p>
836      *  {@code acos(z) = -i (log(z + i (sqrt(1 - z<sup>2</sup>))))}
837      * </p>
838      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
839      * input argument is {@code NaN} or infinite.
840      *
841      * @return the inverse cosine of this complex number.
842      */
843     @Override
844     public FieldComplex<T> acos() {
845         if (isNaN) {
846             return getNaN(getPartsField());
847         }
848 
849         return this.add(this.sqrt1z().multiplyPlusI()).log().multiplyMinusI();
850     }
851 
852     /**
853      * Compute the
854      * <a href="http://mathworld.wolfram.com/InverseSine.html" TARGET="_top">
855      * inverse sine</a> of this complex number.
856      * Implements the formula:
857      * <p>
858      *  {@code asin(z) = -i (log(sqrt(1 - z<sup>2</sup>) + iz))}
859      * </p><p>
860      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
861      * input argument is {@code NaN} or infinite.</p>
862      *
863      * @return the inverse sine of this complex number.
864      */
865     @Override
866     public FieldComplex<T> asin() {
867         if (isNaN) {
868             return getNaN(getPartsField());
869         }
870 
871         return sqrt1z().add(this.multiplyPlusI()).log().multiplyMinusI();
872     }
873 
874     /**
875      * Compute the
876      * <a href="http://mathworld.wolfram.com/InverseTangent.html" TARGET="_top">
877      * inverse tangent</a> of this complex number.
878      * Implements the formula:
879      * <p>
880      * {@code atan(z) = (i/2) log((1 - iz)/(1 + iz))}
881      * </p><p>
882      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
883      * input argument is {@code NaN} or infinite.</p>
884      *
885      * @return the inverse tangent of this complex number
886      */
887     @Override
888     public FieldComplex<T> atan() {
889         if (isNaN) {
890             return getNaN(getPartsField());
891         }
892 
893         final T one = getPartsField().getOne();
894         if (real.isZero()) {
895 
896             // singularity at ±i
897             if (imaginary.multiply(imaginary).subtract(one).isZero()) {
898                 return getNaN(getPartsField());
899             }
900 
901             // branch cut on imaginary axis
902             final T zero = getPartsField().getZero();
903             final FieldComplex<T> tmp = createComplex(one.add(imaginary).divide(one.subtract(imaginary)), zero).
904                                         log().multiplyPlusI().multiply(0.5);
905             return createComplex(FastMath.copySign(tmp.real, real), tmp.imaginary);
906 
907         } else {
908             // regular formula
909             final FieldComplex<T> n = createComplex(one.add(imaginary), real.negate());
910             final FieldComplex<T> d = createComplex(one.subtract(imaginary),  real);
911             return n.divide(d).log().multiplyPlusI().multiply(0.5);
912         }
913 
914     }
915 
916     /**
917      * Compute the
918      * <a href="http://mathworld.wolfram.com/Cosine.html" TARGET="_top">
919      * cosine</a> of this complex number.
920      * Implements the formula:
921      * <p>
922      *  {@code cos(a + bi) = cos(a)cosh(b) - sin(a)sinh(b)i}
923      * </p><p>
924      * where the (real) functions on the right-hand side are
925      * {@link FastMath#sin}, {@link FastMath#cos},
926      * {@link FastMath#cosh} and {@link FastMath#sinh}.
927      * </p><p>
928      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
929      * input argument is {@code NaN}.
930      * </p><p>
931      * Infinite values in real or imaginary parts of the input may result in
932      * infinite or NaN values returned in parts of the result.</p>
933      * <pre>
934      *  Examples:
935      *  <code>
936      *   cos(1 &plusmn; INFINITY i) = 1 \u2213 INFINITY i
937      *   cos(&plusmn;INFINITY + i) = NaN + NaN i
938      *   cos(&plusmn;INFINITY &plusmn; INFINITY i) = NaN + NaN i
939      *  </code>
940      * </pre>
941      *
942      * @return the cosine of this complex number.
943      */
944     @Override
945     public FieldComplex<T> cos() {
946         if (isNaN) {
947             return getNaN(getPartsField());
948         }
949 
950         final FieldSinCos<T>   scr  = FastMath.sinCos(real);
951         final FieldSinhCosh<T> schi = FastMath.sinhCosh(imaginary);
952         return createComplex(scr.cos().multiply(schi.cosh()), scr.sin().negate().multiply(schi.sinh()));
953     }
954 
955     /**
956      * Compute the
957      * <a href="http://mathworld.wolfram.com/HyperbolicCosine.html" TARGET="_top">
958      * hyperbolic cosine</a> of this complex number.
959      * Implements the formula:
960      * <pre>
961      *  <code>
962      *   cosh(a + bi) = cosh(a)cos(b) + sinh(a)sin(b)i
963      *  </code>
964      * </pre>
965      * where the (real) functions on the right-hand side are
966      * {@link FastMath#sin}, {@link FastMath#cos},
967      * {@link FastMath#cosh} and {@link FastMath#sinh}.
968      * <p>
969      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
970      * input argument is {@code NaN}.
971      * </p>
972      * Infinite values in real or imaginary parts of the input may result in
973      * infinite or NaN values returned in parts of the result.
974      * <pre>
975      *  Examples:
976      *  <code>
977      *   cosh(1 &plusmn; INFINITY i) = NaN + NaN i
978      *   cosh(&plusmn;INFINITY + i) = INFINITY &plusmn; INFINITY i
979      *   cosh(&plusmn;INFINITY &plusmn; INFINITY i) = NaN + NaN i
980      *  </code>
981      * </pre>
982      *
983      * @return the hyperbolic cosine of this complex number.
984      */
985     @Override
986     public FieldComplex<T> cosh() {
987         if (isNaN) {
988             return getNaN(getPartsField());
989         }
990 
991         final FieldSinhCosh<T> schr = FastMath.sinhCosh(real);
992         final FieldSinCos<T>   sci  = FastMath.sinCos(imaginary);
993         return createComplex(schr.cosh().multiply(sci.cos()), schr.sinh().multiply(sci.sin()));
994     }
995 
996     /**
997      * Compute the
998      * <a href="http://mathworld.wolfram.com/ExponentialFunction.html" TARGET="_top">
999      * exponential function</a> of this complex number.
1000      * Implements the formula:
1001      * <pre>
1002      *  <code>
1003      *   exp(a + bi) = exp(a)cos(b) + exp(a)sin(b)i
1004      *  </code>
1005      * </pre>
1006      * where the (real) functions on the right-hand side are
1007      * {@link FastMath#exp(CalculusFieldElement)} p}, {@link FastMath#cos(CalculusFieldElement)}, and
1008      * {@link FastMath#sin(CalculusFieldElement)}.
1009      * <p>
1010      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
1011      * input argument is {@code NaN}.
1012      * </p>
1013      * Infinite values in real or imaginary parts of the input may result in
1014      * infinite or NaN values returned in parts of the result.
1015      * <pre>
1016      *  Examples:
1017      *  <code>
1018      *   exp(1 &plusmn; INFINITY i) = NaN + NaN i
1019      *   exp(INFINITY + i) = INFINITY + INFINITY i
1020      *   exp(-INFINITY + i) = 0 + 0i
1021      *   exp(&plusmn;INFINITY &plusmn; INFINITY i) = NaN + NaN i
1022      *  </code>
1023      * </pre>
1024      *
1025      * @return <code><i>e</i><sup>this</sup></code>.
1026      */
1027     @Override
1028     public FieldComplex<T> exp() {
1029         if (isNaN) {
1030             return getNaN(getPartsField());
1031         }
1032 
1033         final T              expReal = FastMath.exp(real);
1034         final FieldSinCos<T> sc      = FastMath.sinCos(imaginary);
1035         return createComplex(expReal.multiply(sc.cos()), expReal.multiply(sc.sin()));
1036     }
1037 
1038     /** {@inheritDoc} */
1039     @Override
1040     public FieldComplex<T> expm1() {
1041         if (isNaN) {
1042             return getNaN(getPartsField());
1043         }
1044 
1045         final T              expm1Real = FastMath.expm1(real);
1046         final FieldSinCos<T> sc        = FastMath.sinCos(imaginary);
1047         return createComplex(expm1Real.multiply(sc.cos()), expm1Real.multiply(sc.sin()));
1048     }
1049 
1050     /**
1051      * Compute the
1052      * <a href="http://mathworld.wolfram.com/NaturalLogarithm.html" TARGET="_top">
1053      * natural logarithm</a> of this complex number.
1054      * Implements the formula:
1055      * <pre>
1056      *  <code>
1057      *   log(a + bi) = ln(|a + bi|) + arg(a + bi)i
1058      *  </code>
1059      * </pre>
1060      * where ln on the right hand side is {@link FastMath#log(CalculusFieldElement)},
1061      * {@code |a + bi|} is the modulus, {@link #abs},  and
1062      * {@code arg(a + bi) = }{@link FastMath#atan2}(b, a).
1063      * <p>
1064      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
1065      * input argument is {@code NaN}.
1066      * </p>
1067      * Infinite (or critical) values in real or imaginary parts of the input may
1068      * result in infinite or NaN values returned in parts of the result.
1069      * <pre>
1070      *  Examples:
1071      *  <code>
1072      *   log(1 &plusmn; INFINITY i) = INFINITY &plusmn; (&pi;/2)i
1073      *   log(INFINITY + i) = INFINITY + 0i
1074      *   log(-INFINITY + i) = INFINITY + &pi;i
1075      *   log(INFINITY &plusmn; INFINITY i) = INFINITY &plusmn; (&pi;/4)i
1076      *   log(-INFINITY &plusmn; INFINITY i) = INFINITY &plusmn; (3&pi;/4)i
1077      *   log(0 + 0i) = -INFINITY + 0i
1078      *  </code>
1079      * </pre>
1080      *
1081      * @return the value <code>ln &nbsp; this</code>, the natural logarithm
1082      * of {@code this}.
1083      */
1084     @Override
1085     public FieldComplex<T> log() {
1086         if (isNaN) {
1087             return getNaN(getPartsField());
1088         }
1089 
1090         return createComplex(FastMath.log(FastMath.hypot(real, imaginary)),
1091                              FastMath.atan2(imaginary, real));
1092     }
1093 
1094     /** {@inheritDoc} */
1095     @Override
1096     public FieldComplex<T> log1p() {
1097         return add(1.0).log();
1098     }
1099 
1100     /** {@inheritDoc} */
1101     @Override
1102     public FieldComplex<T> log10() {
1103         return log().divide(LOG10);
1104     }
1105 
1106     /**
1107      * Returns of value of this complex number raised to the power of {@code x}.
1108      * <p>
1109      * If {@code x} is a real number whose real part has an integer value, returns {@link #pow(int)},
1110      * if both {@code this} and {@code x} are real and {@link FastMath#pow(double, double)}
1111      * with the corresponding real arguments would return a finite number (neither NaN
1112      * nor infinite), then returns the same value converted to {@code Complex},
1113      * with the same special cases.
1114      * In all other cases real cases, implements y<sup>x</sup> = exp(x&middot;log(y)).
1115      * </p>
1116      *
1117      * @param  x exponent to which this {@code Complex} is to be raised.
1118      * @return <code> this<sup>x</sup></code>.
1119      * @throws NullArgumentException if x is {@code null}.
1120      */
1121     @Override
1122     public FieldComplex<T> pow(FieldComplex<T> x)
1123         throws NullArgumentException {
1124 
1125         MathUtils.checkNotNull(x);
1126 
1127         if (x.imaginary.isZero()) {
1128             final int nx = (int) FastMath.rint(x.real.getReal());
1129             if (x.real.getReal() == nx) {
1130                 // integer power
1131                 return pow(nx);
1132             } else if (this.imaginary.isZero()) {
1133                 // check real implementation that handles a bunch of special cases
1134                 final T realPow = FastMath.pow(this.real, x.real);
1135                 if (realPow.isFinite()) {
1136                     return createComplex(realPow, getPartsField().getZero());
1137                 }
1138             }
1139         }
1140 
1141         // generic implementation
1142         return this.log().multiply(x).exp();
1143 
1144     }
1145 
1146 
1147     /**
1148      * Returns of value of this complex number raised to the power of {@code x}.
1149      * <p>
1150      * If {@code x} has an integer value, returns {@link #pow(int)},
1151      * if {@code this} is real and {@link FastMath#pow(double, double)}
1152      * with the corresponding real arguments would return a finite number (neither NaN
1153      * nor infinite), then returns the same value converted to {@code Complex},
1154      * with the same special cases.
1155      * In all other cases real cases, implements y<sup>x</sup> = exp(x&middot;log(y)).
1156      * </p>
1157      *
1158      * @param  x exponent to which this {@code Complex} is to be raised.
1159      * @return <code> this<sup>x</sup></code>.
1160      */
1161     public FieldComplex<T> pow(T x) {
1162 
1163         final int nx = (int) FastMath.rint(x.getReal());
1164         if (x.getReal() == nx) {
1165             // integer power
1166             return pow(nx);
1167         } else if (this.imaginary.isZero()) {
1168             // check real implementation that handles a bunch of special cases
1169             final T realPow = FastMath.pow(this.real, x);
1170             if (realPow.isFinite()) {
1171                 return createComplex(realPow, getPartsField().getZero());
1172             }
1173         }
1174 
1175         // generic implementation
1176         return this.log().multiply(x).exp();
1177 
1178     }
1179 
1180     /**
1181      * Returns of value of this complex number raised to the power of {@code x}.
1182      * <p>
1183      * If {@code x} has an integer value, returns {@link #pow(int)},
1184      * if {@code this} is real and {@link FastMath#pow(double, double)}
1185      * with the corresponding real arguments would return a finite number (neither NaN
1186      * nor infinite), then returns the same value converted to {@code Complex},
1187      * with the same special cases.
1188      * In all other cases real cases, implements y<sup>x</sup> = exp(x&middot;log(y)).
1189      * </p>
1190      *
1191      * @param  x exponent to which this {@code Complex} is to be raised.
1192      * @return <code> this<sup>x</sup></code>.
1193      */
1194     @Override
1195     public FieldComplex<T> pow(double x) {
1196 
1197         final int nx = (int) FastMath.rint(x);
1198         if (x == nx) {
1199             // integer power
1200             return pow(nx);
1201         } else if (this.imaginary.isZero()) {
1202             // check real implementation that handles a bunch of special cases
1203             final T realPow = FastMath.pow(this.real, x);
1204             if (realPow.isFinite()) {
1205                 return createComplex(realPow, getPartsField().getZero());
1206             }
1207         }
1208 
1209         // generic implementation
1210         return this.log().multiply(x).exp();
1211 
1212     }
1213 
1214      /** {@inheritDoc} */
1215     @Override
1216     public FieldComplex<T> pow(final int n) {
1217 
1218         FieldComplex<T> result = getField().getOne();
1219         final boolean invert;
1220         int p = n;
1221         if (p < 0) {
1222             invert = true;
1223             p = -p;
1224         } else {
1225             invert = false;
1226         }
1227 
1228         // Exponentiate by successive squaring
1229         FieldComplex<T> square = this;
1230         while (p > 0) {
1231             if ((p & 0x1) > 0) {
1232                 result = result.multiply(square);
1233             }
1234             square = square.multiply(square);
1235             p = p >> 1;
1236         }
1237 
1238         return invert ? result.reciprocal() : result;
1239 
1240     }
1241 
1242      /**
1243       * Compute the
1244      * <a href="http://mathworld.wolfram.com/Sine.html" TARGET="_top">
1245      * sine</a>
1246      * of this complex number.
1247      * Implements the formula:
1248      * <pre>
1249      *  <code>
1250      *   sin(a + bi) = sin(a)cosh(b) + cos(a)sinh(b)i
1251      *  </code>
1252      * </pre>
1253      * where the (real) functions on the right-hand side are
1254      * {@link FastMath#sin}, {@link FastMath#cos},
1255      * {@link FastMath#cosh} and {@link FastMath#sinh}.
1256      * <p>
1257      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
1258      * input argument is {@code NaN}.
1259      * </p><p>
1260      * Infinite values in real or imaginary parts of the input may result in
1261      * infinite or {@code NaN} values returned in parts of the result.
1262      * <pre>
1263      *  Examples:
1264      *  <code>
1265      *   sin(1 &plusmn; INFINITY i) = 1 &plusmn; INFINITY i
1266      *   sin(&plusmn;INFINITY + i) = NaN + NaN i
1267      *   sin(&plusmn;INFINITY &plusmn; INFINITY i) = NaN + NaN i
1268      *  </code>
1269      * </pre>
1270      *
1271      * @return the sine of this complex number.
1272      */
1273     @Override
1274     public FieldComplex<T> sin() {
1275         if (isNaN) {
1276             return getNaN(getPartsField());
1277         }
1278 
1279         final FieldSinCos<T>   scr  = FastMath.sinCos(real);
1280         final FieldSinhCosh<T> schi = FastMath.sinhCosh(imaginary);
1281         return createComplex(scr.sin().multiply(schi.cosh()), scr.cos().multiply(schi.sinh()));
1282 
1283     }
1284 
1285     /** {@inheritDoc}
1286      */
1287     @Override
1288     public FieldSinCos<FieldComplex<T>> sinCos() {
1289         if (isNaN) {
1290             return new FieldSinCos<>(getNaN(getPartsField()), getNaN(getPartsField()));
1291         }
1292 
1293         final FieldSinCos<T>   scr = FastMath.sinCos(real);
1294         final FieldSinhCosh<T> schi = FastMath.sinhCosh(imaginary);
1295         return new FieldSinCos<>(createComplex(scr.sin().multiply(schi.cosh()), scr.cos().multiply(schi.sinh())),
1296                                  createComplex(scr.cos().multiply(schi.cosh()), scr.sin().negate().multiply(schi.sinh())));
1297     }
1298 
1299     /** {@inheritDoc} */
1300     @Override
1301     public FieldComplex<T> atan2(FieldComplex<T> x) {
1302 
1303         // compute r = sqrt(x^2+y^2)
1304         final FieldComplex<T> r = x.square().add(multiply(this)).sqrt();
1305 
1306         if (x.real.getReal() >= 0) {
1307             // compute atan2(y, x) = 2 atan(y / (r + x))
1308             return divide(r.add(x)).atan().multiply(2);
1309         } else {
1310             // compute atan2(y, x) = +/- pi - 2 atan(y / (r - x))
1311             return divide(r.subtract(x)).atan().multiply(-2).add(x.real.getPi());
1312         }
1313     }
1314 
1315     /** {@inheritDoc}
1316      * <p>
1317      * Branch cuts are on the real axis, below +1.
1318      * </p>
1319      */
1320     @Override
1321     public FieldComplex<T> acosh() {
1322         final FieldComplex<T> sqrtPlus  = add(1).sqrt();
1323         final FieldComplex<T> sqrtMinus = subtract(1).sqrt();
1324         return add(sqrtPlus.multiply(sqrtMinus)).log();
1325     }
1326 
1327     /** {@inheritDoc}
1328      * <p>
1329      * Branch cuts are on the imaginary axis, above +i and below -i.
1330      * </p>
1331      */
1332     @Override
1333     public FieldComplex<T> asinh() {
1334         return add(multiply(this).add(1.0).sqrt()).log();
1335     }
1336 
1337     /** {@inheritDoc}
1338      * <p>
1339      * Branch cuts are on the real axis, above +1 and below -1.
1340      * </p>
1341      */
1342     @Override
1343     public FieldComplex<T> atanh() {
1344         final FieldComplex<T> logPlus  = add(1).log();
1345         final FieldComplex<T> logMinus = createComplex(getPartsField().getOne().subtract(real), imaginary.negate()).log();
1346         return logPlus.subtract(logMinus).multiply(0.5);
1347     }
1348 
1349     /**
1350      * Compute the
1351      * <a href="http://mathworld.wolfram.com/HyperbolicSine.html" TARGET="_top">
1352      * hyperbolic sine</a> of this complex number.
1353      * Implements the formula:
1354      * <pre>
1355      *  <code>
1356      *   sinh(a + bi) = sinh(a)cos(b)) + cosh(a)sin(b)i
1357      *  </code>
1358      * </pre>
1359      * where the (real) functions on the right-hand side are
1360      * {@link FastMath#sin}, {@link FastMath#cos},
1361      * {@link FastMath#cosh} and {@link FastMath#sinh}.
1362      * <p>
1363      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
1364      * input argument is {@code NaN}.
1365      * </p><p>
1366      * Infinite values in real or imaginary parts of the input may result in
1367      * infinite or NaN values returned in parts of the result.
1368      * <pre>
1369      *  Examples:
1370      *  <code>
1371      *   sinh(1 &plusmn; INFINITY i) = NaN + NaN i
1372      *   sinh(&plusmn;INFINITY + i) = &plusmn; INFINITY + INFINITY i
1373      *   sinh(&plusmn;INFINITY &plusmn; INFINITY i) = NaN + NaN i
1374      *  </code>
1375      * </pre>
1376      *
1377      * @return the hyperbolic sine of {@code this}.
1378      */
1379     @Override
1380     public FieldComplex<T> sinh() {
1381         if (isNaN) {
1382             return getNaN(getPartsField());
1383         }
1384 
1385         final FieldSinhCosh<T> schr = FastMath.sinhCosh(real);
1386         final FieldSinCos<T>   sci  = FastMath.sinCos(imaginary);
1387         return createComplex(schr.sinh().multiply(sci.cos()), schr.cosh().multiply(sci.sin()));
1388     }
1389 
1390     /** {@inheritDoc}
1391      */
1392     @Override
1393     public FieldSinhCosh<FieldComplex<T>> sinhCosh() {
1394         if (isNaN) {
1395             return new FieldSinhCosh<>(getNaN(getPartsField()), getNaN(getPartsField()));
1396         }
1397 
1398         final FieldSinhCosh<T> schr = FastMath.sinhCosh(real);
1399         final FieldSinCos<T>   sci  = FastMath.sinCos(imaginary);
1400         return new FieldSinhCosh<>(createComplex(schr.sinh().multiply(sci.cos()), schr.cosh().multiply(sci.sin())),
1401                                    createComplex(schr.cosh().multiply(sci.cos()), schr.sinh().multiply(sci.sin())));
1402     }
1403 
1404     /**
1405      * Compute the
1406      * <a href="http://mathworld.wolfram.com/SquareRoot.html" TARGET="_top">
1407      * square root</a> of this complex number.
1408      * Implements the following algorithm to compute {@code sqrt(a + bi)}:
1409      * <ol><li>Let {@code t = sqrt((|a| + |a + bi|) / 2)}</li>
1410      * <li><pre>if {@code  a ≥ 0} return {@code t + (b/2t)i}
1411      *  else return {@code |b|/2t + sign(b)t i }</pre></li>
1412      * </ol>
1413      * where <ul>
1414      * <li>{@code |a| = }{@link FastMath#abs(CalculusFieldElement) abs(a)}</li>
1415      * <li>{@code |a + bi| = }{@link FastMath#hypot(CalculusFieldElement, CalculusFieldElement) hypot(a, b)}</li>
1416      * <li>{@code sign(b) = }{@link FastMath#copySign(CalculusFieldElement, CalculusFieldElement) copySign(1, b)}
1417      * </ul>
1418      * The real part is therefore always nonnegative.
1419      * <p>
1420      * Returns {@link #getNaN(Field) NaN} if either real or imaginary part of the
1421      * input argument is {@code NaN}.
1422      * </p>
1423      * <p>
1424      * Infinite values in real or imaginary parts of the input may result in
1425      * infinite or NaN values returned in parts of the result.
1426      * </p>
1427      * <pre>
1428      *  Examples:
1429      *  <code>
1430      *   sqrt(1 ± ∞ i) = ∞ + NaN i
1431      *   sqrt(∞ + i) = ∞ + 0i
1432      *   sqrt(-∞ + i) = 0 + ∞ i
1433      *   sqrt(∞ ± ∞ i) = ∞ + NaN i
1434      *   sqrt(-∞ ± ∞ i) = NaN ± ∞ i
1435      *  </code>
1436      * </pre>
1437      *
1438      * @return the square root of {@code this} with nonnegative real part.
1439      */
1440     @Override
1441     public FieldComplex<T> sqrt() {
1442         if (isNaN) {
1443             return getNaN(getPartsField());
1444         }
1445 
1446         if (isZero()) {
1447             return getZero(getPartsField());
1448         }
1449 
1450         T t = FastMath.sqrt((FastMath.abs(real).add(FastMath.hypot(real, imaginary))).multiply(0.5));
1451         if (real.getReal() >= 0.0) {
1452             return createComplex(t, imaginary.divide(t.multiply(2)));
1453         } else {
1454             return createComplex(FastMath.abs(imaginary).divide(t.multiply(2)),
1455                                  FastMath.copySign(t, imaginary));
1456         }
1457     }
1458 
1459     /**
1460      * Compute the
1461      * <a href="http://mathworld.wolfram.com/SquareRoot.html" TARGET="_top">
1462      * square root</a> of <code>1 - this<sup>2</sup></code> for this complex
1463      * number.
1464      * Computes the result directly as
1465      * {@code sqrt(ONE.subtract(z.square()))}.
1466      * <p>
1467      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
1468      * input argument is {@code NaN}.
1469      * </p>
1470      * Infinite values in real or imaginary parts of the input may result in
1471      * infinite or NaN values returned in parts of the result.
1472      *
1473      * @return the square root of <code>1 - this<sup>2</sup></code>.
1474      */
1475     public FieldComplex<T> sqrt1z() {
1476         final FieldComplex<T> t2 = this.square();
1477         return createComplex(getPartsField().getOne().subtract(t2.real), t2.imaginary.negate()).sqrt();
1478     }
1479 
1480     /** {@inheritDoc}
1481      * <p>
1482      * This implementation compute the principal cube root by using a branch cut along real negative axis.
1483      * </p>
1484      */
1485     @Override
1486     public FieldComplex<T> cbrt() {
1487         final T              magnitude = FastMath.cbrt(abs().getRealPart());
1488         final FieldSinCos<T> sc        = FastMath.sinCos(getArgument().divide(3));
1489         return createComplex(magnitude.multiply(sc.cos()), magnitude.multiply(sc.sin()));
1490     }
1491 
1492     /** {@inheritDoc}
1493      * <p>
1494      * This implementation compute the principal n<sup>th</sup> root by using a branch cut along real negative axis.
1495      * </p>
1496      */
1497     @Override
1498     public FieldComplex<T> rootN(int n) {
1499         final T              magnitude = FastMath.pow(abs().getRealPart(), 1.0 / n);
1500         final FieldSinCos<T> sc        = FastMath.sinCos(getArgument().divide(n));
1501         return createComplex(magnitude.multiply(sc.cos()), magnitude.multiply(sc.sin()));
1502     }
1503 
1504     /**
1505      * Compute the
1506      * <a href="http://mathworld.wolfram.com/Tangent.html" TARGET="_top">
1507      * tangent</a> of this complex number.
1508      * Implements the formula:
1509      * <pre>
1510      *  <code>
1511      *   tan(a + bi) = sin(2a)/(cos(2a)+cosh(2b)) + [sinh(2b)/(cos(2a)+cosh(2b))]i
1512      *  </code>
1513      * </pre>
1514      * where the (real) functions on the right-hand side are
1515      * {@link FastMath#sin}, {@link FastMath#cos}, {@link FastMath#cosh} and
1516      * {@link FastMath#sinh}.
1517      * <p>
1518      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
1519      * input argument is {@code NaN}.
1520      * </p>
1521      * Infinite (or critical) values in real or imaginary parts of the input may
1522      * result in infinite or NaN values returned in parts of the result.
1523      * <pre>
1524      *  Examples:
1525      *  <code>
1526      *   tan(a &plusmn; INFINITY i) = 0 &plusmn; i
1527      *   tan(&plusmn;INFINITY + bi) = NaN + NaN i
1528      *   tan(&plusmn;INFINITY &plusmn; INFINITY i) = NaN + NaN i
1529      *   tan(&plusmn;&pi;/2 + 0 i) = &plusmn;INFINITY + NaN i
1530      *  </code>
1531      * </pre>
1532      *
1533      * @return the tangent of {@code this}.
1534      */
1535     @Override
1536     public FieldComplex<T> tan() {
1537         if (isNaN || real.isInfinite()) {
1538             return getNaN(getPartsField());
1539         }
1540         if (imaginary.getReal() > 20.0) {
1541             return getI(getPartsField());
1542         }
1543         if (imaginary.getReal() < -20.0) {
1544             return getMinusI(getPartsField());
1545         }
1546 
1547         final FieldSinCos<T> sc2r = FastMath.sinCos(real.multiply(2));
1548         T imaginary2 = imaginary.multiply(2);
1549         T d = sc2r.cos().add(FastMath.cosh(imaginary2));
1550 
1551         return createComplex(sc2r.sin().divide(d), FastMath.sinh(imaginary2).divide(d));
1552 
1553     }
1554 
1555     /**
1556      * Compute the
1557      * <a href="http://mathworld.wolfram.com/HyperbolicTangent.html" TARGET="_top">
1558      * hyperbolic tangent</a> of this complex number.
1559      * Implements the formula:
1560      * <pre>
1561      *  <code>
1562      *   tan(a + bi) = sinh(2a)/(cosh(2a)+cos(2b)) + [sin(2b)/(cosh(2a)+cos(2b))]i
1563      *  </code>
1564      * </pre>
1565      * where the (real) functions on the right-hand side are
1566      * {@link FastMath#sin}, {@link FastMath#cos}, {@link FastMath#cosh} and
1567      * {@link FastMath#sinh}.
1568      * <p>
1569      * Returns {@link #getNaN(Field)} if either real or imaginary part of the
1570      * input argument is {@code NaN}.
1571      * </p>
1572      * Infinite values in real or imaginary parts of the input may result in
1573      * infinite or NaN values returned in parts of the result.
1574      * <pre>
1575      *  Examples:
1576      *  <code>
1577      *   tanh(a &plusmn; INFINITY i) = NaN + NaN i
1578      *   tanh(&plusmn;INFINITY + bi) = &plusmn;1 + 0 i
1579      *   tanh(&plusmn;INFINITY &plusmn; INFINITY i) = NaN + NaN i
1580      *   tanh(0 + (&pi;/2)i) = NaN + INFINITY i
1581      *  </code>
1582      * </pre>
1583      *
1584      * @return the hyperbolic tangent of {@code this}.
1585      */
1586     @Override
1587     public FieldComplex<T> tanh() {
1588         if (isNaN || imaginary.isInfinite()) {
1589             return getNaN(getPartsField());
1590         }
1591         if (real.getReal() > 20.0) {
1592             return getOne(getPartsField());
1593         }
1594         if (real.getReal() < -20.0) {
1595             return getMinusOne(getPartsField());
1596         }
1597         T real2 = real.multiply(2);
1598         final FieldSinCos<T> sc2i = FastMath.sinCos(imaginary.multiply(2));
1599         T d = FastMath.cosh(real2).add(sc2i.cos());
1600 
1601         return createComplex(FastMath.sinh(real2).divide(d), sc2i.sin().divide(d));
1602     }
1603 
1604 
1605 
1606     /**
1607      * Compute the argument of this complex number.
1608      * The argument is the angle phi between the positive real axis and
1609      * the point representing this number in the complex plane.
1610      * The value returned is between -PI (not inclusive)
1611      * and PI (inclusive), with negative values returned for numbers with
1612      * negative imaginary parts.
1613      * <p>
1614      * If either real or imaginary part (or both) is NaN, NaN is returned.
1615      * Infinite parts are handled as {@code Math.atan2} handles them,
1616      * essentially treating finite parts as zero in the presence of an
1617      * infinite coordinate and returning a multiple of pi/4 depending on
1618      * the signs of the infinite parts.
1619      * See the javadoc for {@code Math.atan2} for full details.
1620      *
1621      * @return the argument of {@code this}.
1622      */
1623     public T getArgument() {
1624         return FastMath.atan2(getImaginaryPart(), getRealPart());
1625     }
1626 
1627     /**
1628      * Computes the n-th roots of this complex number.
1629      * The nth roots are defined by the formula:
1630      * <pre>
1631      *  <code>
1632      *   z<sub>k</sub> = abs<sup>1/n</sup> (cos(phi + 2&pi;k/n) + i (sin(phi + 2&pi;k/n))
1633      *  </code>
1634      * </pre>
1635      * for <i>{@code k=0, 1, ..., n-1}</i>, where {@code abs} and {@code phi}
1636      * are respectively the {@link #abs() modulus} and
1637      * {@link #getArgument() argument} of this complex number.
1638      * <p>
1639      * If one or both parts of this complex number is NaN, a list with just
1640      * one element, {@link #getNaN(Field)} is returned.
1641      * if neither part is NaN, but at least one part is infinite, the result
1642      * is a one-element list containing {@link #getInf(Field)}.
1643      *
1644      * @param n Degree of root.
1645      * @return a List of all {@code n}-th roots of {@code this}.
1646      * @throws MathIllegalArgumentException if {@code n <= 0}.
1647      */
1648     public List<FieldComplex<T>> nthRoot(int n) throws MathIllegalArgumentException {
1649 
1650         if (n <= 0) {
1651             throw new MathIllegalArgumentException(LocalizedCoreFormats.CANNOT_COMPUTE_NTH_ROOT_FOR_NEGATIVE_N,
1652                                                    n);
1653         }
1654 
1655         final List<FieldComplex<T>> result = new ArrayList<>();
1656 
1657         if (isNaN) {
1658             result.add(getNaN(getPartsField()));
1659             return result;
1660         }
1661         if (isInfinite()) {
1662             result.add(getInf(getPartsField()));
1663             return result;
1664         }
1665 
1666         // nth root of abs -- faster / more accurate to use a solver here?
1667         final T nthRootOfAbs = FastMath.pow(FastMath.hypot(real, imaginary), 1.0 / n);
1668 
1669         // Compute nth roots of complex number with k = 0, 1, ... n-1
1670         final T nthPhi = getArgument().divide(n);
1671         final double slice = 2 * FastMath.PI / n;
1672         T innerPart = nthPhi;
1673         for (int k = 0; k < n ; k++) {
1674             // inner part
1675             final FieldSinCos<T> scInner = FastMath.sinCos(innerPart);
1676             final T realPart = nthRootOfAbs.multiply(scInner.cos());
1677             final T imaginaryPart = nthRootOfAbs.multiply(scInner.sin());
1678             result.add(createComplex(realPart, imaginaryPart));
1679             innerPart = innerPart.add(slice);
1680         }
1681 
1682         return result;
1683     }
1684 
1685     /**
1686      * Create a complex number given the real and imaginary parts.
1687      *
1688      * @param realPart Real part.
1689      * @param imaginaryPart Imaginary part.
1690      * @return a new complex number instance.
1691      *
1692      * @see #valueOf(CalculusFieldElement, CalculusFieldElement)
1693      */
1694     protected FieldComplex<T> createComplex(final T realPart, final T imaginaryPart) {
1695         return new FieldComplex<>(realPart, imaginaryPart);
1696     }
1697 
1698     /**
1699      * Create a complex number given the real and imaginary parts.
1700      *
1701      * @param realPart Real part.
1702      * @param imaginaryPart Imaginary part.
1703      * @param <T> the type of the field elements
1704      * @return a Complex instance.
1705      */
1706     public static <T extends CalculusFieldElement<T>> FieldComplex<T>
1707         valueOf(T realPart, T imaginaryPart) {
1708         if (realPart.isNaN() || imaginaryPart.isNaN()) {
1709             return getNaN(realPart.getField());
1710         }
1711         return new FieldComplex<>(realPart, imaginaryPart);
1712     }
1713 
1714     /**
1715      * Create a complex number given only the real part.
1716      *
1717      * @param realPart Real part.
1718      * @param <T> the type of the field elements
1719      * @return a Complex instance.
1720      */
1721     public static <T extends CalculusFieldElement<T>> FieldComplex<T>
1722         valueOf(T realPart) {
1723         if (realPart.isNaN()) {
1724             return getNaN(realPart.getField());
1725         }
1726         return new FieldComplex<>(realPart);
1727     }
1728 
1729     /** {@inheritDoc} */
1730     @Override
1731     public FieldComplex<T> newInstance(double realPart) {
1732         return valueOf(getPartsField().getZero().newInstance(realPart));
1733     }
1734 
1735     /** {@inheritDoc} */
1736     @Override
1737     public FieldComplexField<T> getField() {
1738         return FieldComplexField.getField(getPartsField());
1739     }
1740 
1741     /** Get the {@link Field} the real and imaginary parts belong to.
1742      * @return {@link Field} the real and imaginary parts belong to
1743      */
1744     public Field<T> getPartsField() {
1745         return real.getField();
1746     }
1747 
1748     /** {@inheritDoc} */
1749     @Override
1750     public String toString() {
1751         return "(" + real + ", " + imaginary + ")";
1752     }
1753 
1754     /** {@inheritDoc} */
1755     @Override
1756     public FieldComplex<T> scalb(int n) {
1757         return createComplex(FastMath.scalb(real, n), FastMath.scalb(imaginary, n));
1758     }
1759 
1760     /** {@inheritDoc} */
1761     @Override
1762     public FieldComplex<T> ulp() {
1763         return createComplex(FastMath.ulp(real), FastMath.ulp(imaginary));
1764     }
1765 
1766     /** {@inheritDoc} */
1767     @Override
1768     public FieldComplex<T> hypot(FieldComplex<T> y) {
1769         if (isInfinite() || y.isInfinite()) {
1770             return getInf(getPartsField());
1771         } else if (isNaN() || y.isNaN()) {
1772             return getNaN(getPartsField());
1773         } else {
1774             return square().add(y.square()).sqrt();
1775         }
1776     }
1777 
1778     /** {@inheritDoc} */
1779     @Override
1780     public FieldComplex<T> linearCombination(final FieldComplex<T>[] a, final FieldComplex<T>[] b)
1781         throws MathIllegalArgumentException {
1782         final int n = 2 * a.length;
1783         final T[] realA      = MathArrays.buildArray(getPartsField(), n);
1784         final T[] realB      = MathArrays.buildArray(getPartsField(), n);
1785         final T[] imaginaryA = MathArrays.buildArray(getPartsField(), n);
1786         final T[] imaginaryB = MathArrays.buildArray(getPartsField(), n);
1787         for (int i = 0; i < a.length; ++i)  {
1788             final FieldComplex<T> ai = a[i];
1789             final FieldComplex<T> bi = b[i];
1790             realA[2 * i    ]      = ai.real;
1791             realA[2 * i + 1]      = ai.imaginary.negate();
1792             realB[2 * i    ]      = bi.real;
1793             realB[2 * i + 1]      = bi.imaginary;
1794             imaginaryA[2 * i    ] = ai.real;
1795             imaginaryA[2 * i + 1] = ai.imaginary;
1796             imaginaryB[2 * i    ] = bi.imaginary;
1797             imaginaryB[2 * i + 1] = bi.real;
1798         }
1799         return createComplex(real.linearCombination(realA,  realB),
1800                              real.linearCombination(imaginaryA, imaginaryB));
1801     }
1802 
1803     /** {@inheritDoc} */
1804     @Override
1805     public FieldComplex<T> linearCombination(final double[] a, final FieldComplex<T>[] b)
1806         throws MathIllegalArgumentException {
1807         final int n = a.length;
1808         final T[] realB      = MathArrays.buildArray(getPartsField(), n);
1809         final T[] imaginaryB = MathArrays.buildArray(getPartsField(), n);
1810         for (int i = 0; i < a.length; ++i)  {
1811             final FieldComplex<T> bi = b[i];
1812             realB[i]      = bi.real;
1813             imaginaryB[i] = bi.imaginary;
1814         }
1815         return createComplex(real.linearCombination(a,  realB),
1816                              real.linearCombination(a, imaginaryB));
1817     }
1818 
1819     /** {@inheritDoc} */
1820     @Override
1821     public FieldComplex<T> linearCombination(final FieldComplex<T> a1, final FieldComplex<T> b1, final FieldComplex<T> a2, final FieldComplex<T> b2) {
1822         return createComplex(real.linearCombination(a1.real, b1.real,
1823                                                     a1.imaginary.negate(), b1.imaginary,
1824                                                     a2.real, b2.real,
1825                                                     a2.imaginary.negate(), b2.imaginary),
1826                              real.linearCombination(a1.real, b1.imaginary,
1827                                                     a1.imaginary, b1.real,
1828                                                     a2.real, b2.imaginary,
1829                                                     a2.imaginary, b2.real));
1830     }
1831 
1832     /** {@inheritDoc} */
1833     @Override
1834     public FieldComplex<T> linearCombination(final double a1, final FieldComplex<T> b1, final double a2, final FieldComplex<T> b2) {
1835         return createComplex(real.linearCombination(a1, b1.real,
1836                                                     a2, b2.real),
1837                              real.linearCombination(a1, b1.imaginary,
1838                                                     a2, b2.imaginary));
1839     }
1840 
1841     /** {@inheritDoc} */
1842     @Override
1843     public FieldComplex<T> linearCombination(final FieldComplex<T> a1, final FieldComplex<T> b1,
1844                                                 final FieldComplex<T> a2, final FieldComplex<T> b2,
1845                                                 final FieldComplex<T> a3, final FieldComplex<T> b3) {
1846         FieldComplex<T>[] a = MathArrays.buildArray(getField(), 3);
1847         a[0] = a1;
1848         a[1] = a2;
1849         a[2] = a3;
1850         FieldComplex<T>[] b = MathArrays.buildArray(getField(), 3);
1851         b[0] = b1;
1852         b[1] = b2;
1853         b[2] = b3;
1854         return linearCombination(a, b);
1855     }
1856 
1857     /** {@inheritDoc} */
1858     @Override
1859     public FieldComplex<T> linearCombination(final double a1, final FieldComplex<T> b1,
1860                                                 final double a2, final FieldComplex<T> b2,
1861                                                 final double a3, final FieldComplex<T> b3) {
1862         FieldComplex<T>[] b = MathArrays.buildArray(getField(), 3);
1863         b[0] = b1;
1864         b[1] = b2;
1865         b[2] = b3;
1866         return linearCombination(new double[]  { a1, a2, a3 }, b);
1867     }
1868 
1869     /** {@inheritDoc} */
1870     @Override
1871     public FieldComplex<T> linearCombination(final FieldComplex<T> a1, final FieldComplex<T> b1,
1872                                                 final FieldComplex<T> a2, final FieldComplex<T> b2,
1873                                                 final FieldComplex<T> a3, final FieldComplex<T> b3,
1874                                                 final FieldComplex<T> a4, final FieldComplex<T> b4) {
1875         FieldComplex<T>[] a = MathArrays.buildArray(getField(), 4);
1876         a[0] = a1;
1877         a[1] = a2;
1878         a[2] = a3;
1879         a[3] = a4;
1880         FieldComplex<T>[] b = MathArrays.buildArray(getField(), 4);
1881         b[0] = b1;
1882         b[1] = b2;
1883         b[2] = b3;
1884         b[3] = b4;
1885         return linearCombination(a, b);
1886     }
1887 
1888     /** {@inheritDoc} */
1889     @Override
1890     public FieldComplex<T> linearCombination(final double a1, final FieldComplex<T> b1,
1891                                                 final double a2, final FieldComplex<T> b2,
1892                                                 final double a3, final FieldComplex<T> b3,
1893                                                 final double a4, final FieldComplex<T> b4) {
1894         FieldComplex<T>[] b = MathArrays.buildArray(getField(), 4);
1895         b[0] = b1;
1896         b[1] = b2;
1897         b[2] = b3;
1898         b[3] = b4;
1899         return linearCombination(new double[]  { a1, a2, a3, a4 }, b);
1900     }
1901 
1902     /** {@inheritDoc} */
1903     @Override
1904     public FieldComplex<T> ceil() {
1905         return createComplex(FastMath.ceil(getRealPart()), FastMath.ceil(getImaginaryPart()));
1906     }
1907 
1908     /** {@inheritDoc} */
1909     @Override
1910     public FieldComplex<T> floor() {
1911         return createComplex(FastMath.floor(getRealPart()), FastMath.floor(getImaginaryPart()));
1912     }
1913 
1914     /** {@inheritDoc} */
1915     @Override
1916     public FieldComplex<T> rint() {
1917         return createComplex(FastMath.rint(getRealPart()), FastMath.rint(getImaginaryPart()));
1918     }
1919 
1920     /** {@inheritDoc}
1921      * <p>
1922      * for complex numbers, the integer n corresponding to {@code this.subtract(remainder(a)).divide(a)}
1923      * is a <a href="https://en.wikipedia.org/wiki/Gaussian_integer">Wikipedia - Gaussian integer</a>.
1924      * </p>
1925      */
1926     @Override
1927     public FieldComplex<T> remainder(final double a) {
1928         return createComplex(FastMath.IEEEremainder(getRealPart(), a), FastMath.IEEEremainder(getImaginaryPart(), a));
1929     }
1930 
1931     /** {@inheritDoc}
1932      * <p>
1933      * for complex numbers, the integer n corresponding to {@code this.subtract(remainder(a)).divide(a)}
1934      * is a <a href="https://en.wikipedia.org/wiki/Gaussian_integer">Wikipedia - Gaussian integer</a>.
1935      * </p>
1936      */
1937     @Override
1938     public FieldComplex<T> remainder(final FieldComplex<T> a) {
1939         final FieldComplex<T> complexQuotient = divide(a);
1940         final T  qRInt           = FastMath.rint(complexQuotient.real);
1941         final T  qIInt           = FastMath.rint(complexQuotient.imaginary);
1942         return createComplex(real.subtract(qRInt.multiply(a.real)).add(qIInt.multiply(a.imaginary)),
1943                              imaginary.subtract(qRInt.multiply(a.imaginary)).subtract(qIInt.multiply(a.real)));
1944     }
1945 
1946     /** {@inheritDoc} */
1947     @Override
1948     public FieldComplex<T> sign() {
1949         if (isNaN() || isZero()) {
1950             return this;
1951         } else {
1952             return this.divide(FastMath.hypot(real, imaginary));
1953         }
1954     }
1955 
1956     /** {@inheritDoc}
1957      * <p>
1958      * The signs of real and imaginary parts are copied independently.
1959      * </p>
1960      */
1961     @Override
1962     public FieldComplex<T> copySign(final FieldComplex<T> z) {
1963         return createComplex(FastMath.copySign(getRealPart(), z.getRealPart()),
1964                              FastMath.copySign(getImaginaryPart(), z.getImaginaryPart()));
1965     }
1966 
1967     /** {@inheritDoc} */
1968     @Override
1969     public FieldComplex<T> copySign(double r) {
1970         return createComplex(FastMath.copySign(getRealPart(), r), FastMath.copySign(getImaginaryPart(), r));
1971     }
1972 
1973     /** {@inheritDoc} */
1974     @Override
1975     public FieldComplex<T> toDegrees() {
1976         return createComplex(FastMath.toDegrees(getRealPart()), FastMath.toDegrees(getImaginaryPart()));
1977     }
1978 
1979     /** {@inheritDoc} */
1980     @Override
1981     public FieldComplex<T> toRadians() {
1982         return createComplex(FastMath.toRadians(getRealPart()), FastMath.toRadians(getImaginaryPart()));
1983     }
1984 
1985     /** {@inheritDoc} */
1986     @Override
1987     public FieldComplex<T> getPi() {
1988         return getPi(getPartsField());
1989     }
1990 
1991 }