HILA
Loading...
Searching...
No Matches
array.h
Go to the documentation of this file.
1/**
2 * @file array.h
3 * @brief Definition of Array class
4 * @details This file contains the definitions of Array class and utility functions related to it.
5 */
6
7#ifndef ARRAY_H_
8#define ARRAY_H_
9
10#include "datatypes/matrix.h"
11
12/**
13 * @brief \f$ n\times m \f$ Array type
14 * @details Acts as array class which stores data in a simple C style array.
15 *
16 * Main functionality which the Array class offers is to supplement the fall backs of storing
17 * information in a Matrix data structure.
18 *
19 * For example assigning a value to each element with the Matrix class is not directly possible
20 * using the assignment operator=. This is because assignment with matrices is defined as
21 * \f$ M = a = a*I \f$ where M is a general matrix, a is a scalar and I an identity matrix. The
22 * result would only assign a to the diagonal elements.
23 *
24 * Unlike the Matrix object, the Array object assigns element wise. This allows filling the Array
25 * with the assignment operator. Additionally element wise operations are useable as long as they
26 * are defined for the Array type. For example the operation:
27 *
28 * \code
29 * Array<n,m,double> m.random();
30 * sin(m);
31 * \endcode
32 *
33 * is defined, since the \f$sin\f$ function is defined for doubles.
34 *
35 * The above operation would not work for matrices, but with casting operations. Matrix::asArray and
36 * Array::asMatrix one can take full advantage of element wise operations.
37 *
38 * @tparam n Number of rows
39 * @tparam m Number of columns
40 * @tparam T Array element type
41 */
42template <const int n, const int m, typename T>
43class Array {
44 public:
46 "Array requires Complex or arithmetic type");
47
48 // Store Array contents in one dimensional C style array
49 T c[n * m];
50
51 // std incantation for field types
52 using base_type = hila::arithmetic_type<T>;
53 using argument_type = T;
54
55 /**
56 * @brief Default constructor
57 * @details The following ways of constructing a matrix are:
58 *
59 * Allocates undefined \f$ n\times m\f$ Array.
60 *
61 * \code{.cpp}
62 * Array<n,m,MyType> A;
63 * \endcode
64 *
65 */
66 Array() = default;
67 ~Array() = default;
68 /**
69 * @brief Copy constructor
70 * @details
71 * Construction from previously defined Array.
72 *
73 * \code{.cpp}
74 * Array<n,m,MyOtherType> A_0;
75 * //
76 * // A_0 is filled with content
77 * //
78 * Array<n,m,MyType> A = A_0;
79 * \endcode
80 *
81 * @param v Array to assign
82 */
83 Array(const Array<n, m, T> &v) = default;
84
85 /**
86 * @brief Scalar constructor
87 * @details
88 *
89 * Construct with given scalar which is assigned to all elements in array
90 *
91 * \code{.cpp}
92 * MyType x = hila::random();
93 * Array<n,m,MyType> A = x;
94 * \endcode
95 *
96 * @tparam S Type for scalar
97 * @param rhs Scalar to assign
98 */
99 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value, int> = 0>
100 explicit inline Array(const S rhs) {
101 for (int i = 0; i < n * m; i++) {
102 this->c[i] = rhs;
103 }
104 }
105
106 /**
107 * @brief Zero constructor
108 * @details Constructing from 0 sets the whole Array to zero
109 * @internal
110 * @param z 0
111 */
112#pragma hila loop_function
113 inline Array(const std::nullptr_t &z) {
114 for (int i = 0; i < n * m; i++)
115 c[i] = 0;
116 }
117
118 // Construct array automatically from right-size initializer list
119 // This does not seem to be dangerous, so keep non-explicit
120
121 /**
122 * @brief Initializer list constructor
123 *
124 * Construction from c++ initializer list.
125 *
126 * \code{.cpp}
127 * Array<2,2,int> A = {1, 0
128 * 0, 1};
129 * \endcode
130 * @tparam S Element type of initializer list
131 * @param rhs Initializer list to assign
132 */
133 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value, int> = 0>
134 inline Array(std::initializer_list<S> rhs) {
135 assert(rhs.size() == n * m && "Array initializer list size must match variable size");
136 int i = 0;
137 for (auto it = rhs.begin(); it != rhs.end(); it++, i++) {
138 c[i] = *it;
139 }
140 }
141
142 /**
143 * @brief Returns number of rows
144 *
145 * @return constexpr int
146 */
147 constexpr int rows() const {
148 return n;
149 }
150 /**
151 * @brief Returns number of columns
152 *
153 * @return constexpr int
154 */
155 constexpr int columns() const {
156 return m;
157 }
158
159 /**
160 * @brief Returns size of #Vector Array or square Array
161 *
162 * @tparam q Number of rows n
163 * @tparam p Number of columns m
164 * @return constexpr int
165 */
166 template <int q = n, int p = m, std::enable_if_t<q == 1, int> = 0>
167 constexpr int size() const {
168 return p;
169 }
170
171 template <int q = n, int p = m, std::enable_if_t<p == 1, int> = 0>
172 constexpr int size() const {
173 return q;
174 }
175
176 template <int q = n, int p = m, std::enable_if_t<q == p, int> = 0>
177 constexpr int size() const {
178 return q;
179 }
180
181 /**
182 * @brief Standard array indexing operation for 2D Array
183 *
184 * @details Accessing singular elements is insufficient, but matrix elements are often quite
185 * small.
186 *
187 * Example for 2D Array:
188 * \code
189 * Array<n,m,MyType> A;
190 * MyType a = A.e(i,j); // i <= n, j <= m
191 * \endcode
192 *
193 * @param i row index
194 * @param j column index
195 * @return T matrix element type
196 */
197 inline T e(const int i, const int j) const {
198 return c[i * m + j];
199 }
200
201 /// @brief Const overload
202 /// @internal
203 inline T &e(const int i, const int j) const_function {
204 return c[i * m + j];
205 }
206
207 /// @brief Const overload
208 /// @internal
209 template <int q = n, int p = m, std::enable_if_t<(q == 1 || p == 1), int> = 0>
210 inline T e(const int i) const {
211 return c[i];
212 }
213
214 /**
215 * @brief Standard array indexing operation for 1D Array
216 * @details Example for Vector Array:
217 * \code {.cpp}
218 * Array1d<n,MyType> A;
219 * MyType a = A.e(i); // i <= n
220 * \endcode
221 * @param i Vector index type
222 * @return T matrix element type
223 */
224 template <int q = n, int p = m, std::enable_if_t<(q == 1 || p == 1), int> = 0>
225 inline T &e(const int i) const_function {
226 return c[i];
227 }
228
229 /**
230 * @brief Cast Array to Matrix
231 *
232 * @return Matrix<n, m, T>&
233 */
234 Matrix<n, m, T> &asMatrix() const_function {
235 return *reinterpret_cast<Matrix<n, m, T> *>(this);
236 }
237
238 const Matrix<n, m, T> &asMatrix() const {
239 return *reinterpret_cast<const Matrix<n, m, T> *>(this);
240 }
241
242 /**
243 * @brief Cast Array1D to Vector
244 * @details asMatrix will perform the same operation.
245 *
246 * @return Vector<n, T>&
247 */
248 Vector<n, T> &asVector() const_function {
249 static_assert(1 == m, "asVector() only for column arrays");
250 return *reinterpret_cast<Vector<n, T> *>(this);
251 }
252
253 const Vector<n, T> &asVector() const {
254 static_assert(1 == m, "asVector() only for column arrays");
255 return *reinterpret_cast<const Vector<n, T> *>(this);
256 }
257
258 DiagonalMatrix<n, T> &asDiagonalMatrix() const_function {
259 static_assert(1 == m, "asDiagonalMatrix() only for column arrays");
260 return *reinterpret_cast<DiagonalMatrix<n, T> *>(this);
261 }
262
263 const DiagonalMatrix<n, T> &asDiagonalMatrix() const {
264 static_assert(1 == m, "asDiagonalMatrix() only for column arrays");
265 return *reinterpret_cast<const DiagonalMatrix<n, T> *>(this);
266 }
267
268
269 // casting from one Array (number) type to another
270 // TODO: CHECK AVX CONVERSIONS
271 template <typename S, std::enable_if_t<hila::is_assignable<S &, T>::value, int> = 0>
272 operator Array<n, m, S>() {
273 Array<n, m, S> res;
274 for (int i = 0; i < n * m; i++) {
275 res.c[i] = c[i];
276 }
277 return res;
278 }
279
280 /**
281 * @brief Unary - operator
282 * @details Returns Array with the signs of all the elements in the Arrat flipped.
283 *
284 * @return Array<n,m,T>
285 */
286 inline Array<n, m, T> operator-() const {
287 Array<n, m, T> res;
288 for (int i = 0; i < n * m; i++) {
289 res.c[i] = -c[i];
290 }
291 return res;
292 }
293
294 /**
295 * @brief Unary + operator
296 * @details Equivalent to identity operator meaning that Arrat stays as is.
297 *
298 * @return const Array<n, m, T>
299 */
300 inline Array<n, m, T> operator+() const {
301 return *this;
302 }
303
304 /**
305 * @brief Scalar assignment operator
306 * @details The following ways to assign an Array are:
307 *
308 * __Assignment from scalar__:
309 *
310 * Assignment from scalar assigns the scalar uniformly to the array
311 *
312 * \code {.cpp}
313 * MyType a = hila::random;
314 * Array<n,m,MyType> A;
315 * A = a;
316 * \endcode
317 *
318 */
319#pragma hila loop_function
320 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value, int> = 0>
321 inline Array<n, m, T> &operator=(const S rhs) {
322 for (int i = 0; i < n * m; i++) {
323 c[i] = rhs;
324 }
325 return *this;
326 }
327
328 /**
329 * @brief Compare equality of Arrays
330 *
331 * Two Arrays are equal iff Arrays are of same dimension and all elements compare to equal
332 * Note: Complex == scalar if arithmetic value is equal
333 *
334 * @tparam S
335 * @param rhs
336 * @return true if equal
337 */
338
339 template <typename S, int n1, int m1>
340 bool operator==(const Array<n1, m1, S> &rhs) const {
341 if constexpr (n != n1 || m != m1)
342 return false;
343
344 for (int i = 0; i < n; i++)
345 for (int j = 0; j < m; j++) {
346 if (e(i, j) != rhs.e(i, j))
347 return false;
348 }
349 return true;
350 }
351
352 /**
353 * @brief Compare non-equality of two Arrays
354 *
355 * Negation of operator==()
356 *
357 * @tparam S
358 * @tparam n1
359 * @tparam m1
360 * @param rhs
361 * @return true
362 * @return false
363 */
364
365 template <typename S, int n1, int m1>
366 bool operator!=(const Array<n1, m1, S> &rhs) const {
367 return !(*this == rhs);
368 }
369
370 /**
371 * @brief Add assign operator with Array or scalar
372 * @details Add assign operator can be used in the following ways
373 *
374 * __Add assign Array__:
375 *
376 * Requires that Arrays have same dimensions
377 *
378 * \code {.cpp}
379 * Array<n,m,MyType> A,B;
380 * A = 1;
381 * B = 1;
382 * A += B; \\ A is uniformly 2
383 * \endcode
384 *
385 * __Add assign scalar__:
386 *
387 * Adds scalar \f$ a \f$ to Array uniformly
388 *
389 * \code {.cpp}
390 * Array<n,m,MyType> A = 1;
391 * A += 1 ; \\ A is uniformly 2
392 * \endcode
393 *
394 * @tparam S Element type of rhs
395 * @param rhs Array to add
396 * @return Array<n, m, T>
397 */
398#pragma hila loop_function
399 template <typename S, std::enable_if_t<std::is_convertible<S, T>::value, int> = 0>
401 for (int i = 0; i < n * m; i++) {
402 c[i] += rhs.c[i];
403 }
404 return *this;
405 }
406
407 /**
408 * @brief Subtract assign operator with Array or scalar
409 * @details Subtract assign operator can be used in the following ways
410 *
411 * __Subtract assign Array__:
412 *
413 * \code {.cpp}
414 * Array<n,m,MyType> A,B;
415 * A = 3;
416 * B = 1;
417 * A -= B; \\ A is uniformly 2
418 * \endcode
419 *
420 * __Subtract assign scalar__:
421 *
422 * Subtract scalar uniformly from square matrix
423 *
424 * \code {.cpp}
425 * A<n,m,MyType> A = 3;
426 * A -= 1 ; \\ A is uniformly 2
427 * \endcode
428 * @tparam S
429 * @param rhs
430 * @return Array<n, m, T>&
431 */
432 template <typename S, std::enable_if_t<std::is_convertible<S, T>::value, int> = 0>
434 for (int i = 0; i < n * m; i++) {
435 c[i] -= rhs.c[i];
436 }
437 return *this;
438 }
439
440 // add assign type T and convertible
441 template <typename S,
442 std::enable_if_t<std::is_convertible<hila::type_plus<T, S>, T>::value, int> = 0>
443 Array<n, m, T> &operator+=(const S rhs) {
444 for (int i = 0; i < n * m; i++) {
445 c[i] += rhs;
446 }
447 return *this;
448 }
449
450 // subtract assign type T and convertible
451 template <typename S,
452 std::enable_if_t<std::is_convertible<hila::type_minus<T, S>, T>::value, int> = 0>
453 Array<n, m, T> &operator-=(const S rhs) {
454 for (int i = 0; i < n * m; i++) {
455 c[i] -= rhs;
456 }
457 return *this;
458 }
459
460 /**
461 * @brief Multiply assign scalar or array
462 * @details Multiplication works element wise
463 *
464 * Multiply assign operator can be used in the following ways
465 *
466 * __Multiply assign Array__:
467 *
468 * \code {.cpp}
469 * Array<n,m,MyType> A,B;
470 * A = 2;
471 * B = 2;
472 * A *= B; \\ A is uniformly 4
473 * \endcode
474 *
475 * __Multiply assign scalar__:
476 *
477 * \code {.cpp}
478 * Array<n,m,MyType> A = 1;
479
480 * A *= 2 ; \\ A is uniformly 2
481 * \endcode
482 * @tparam S
483 * @param rhs
484 * @return Array<n, m, T>&
485 */
486 template <typename S,
487 std::enable_if_t<std::is_convertible<hila::type_mul<T, S>, T>::value, int> = 0>
489 for (int i = 0; i < n * m; i++) {
490 c[i] *= rhs.c[i];
491 }
492 return *this;
493 }
494
495 /// multiply assign with scalar
496 template <typename S,
497 std::enable_if_t<std::is_convertible<hila::type_mul<T, S>, T>::value, int> = 0>
499 for (int i = 0; i < n * m; i++) {
500 c[i] *= rhs;
501 }
502 return *this;
503 }
504
505 /**
506 * @brief Division assign with array or scalar
507 * @details Division works element wise
508 *
509 * Division assign operator can be used in the following ways
510 *
511 * __Division assign Array__:
512 *
513 * \code {.cpp}
514 * Array<n,m,MyType> A,B;
515 * A = 2;
516 * B = 2;
517 * A /= B; \\ A is uniformly 1
518 * \endcode
519 *
520 * __Division assign scalar__:
521 *
522 * \code {.cpp}
523 * Array<n,m,MyType> A = 2;
524
525 * A /= 2 ; \\ A is uniformly 1
526 * \endcode
527 * @tparam S
528 * @param rhs
529 * @return Array<n, m, T>&
530 */
531 template <typename S,
532 std::enable_if_t<std::is_convertible<hila::type_div<T, S>, T>::value, int> = 0>
534 for (int i = 0; i < n * m; i++) {
535 c[i] /= rhs.c[i];
536 }
537 return *this;
538 }
539
540 // divide assign with scalar
541 template <typename S,
542 std::enable_if_t<std::is_convertible<hila::type_div<T, S>, T>::value, int> = 0>
543 Array<n, m, T> &operator/=(const S rhs) {
544 for (int i = 0; i < n * m; i++) {
545 c[i] /= rhs;
546 }
547 return *this;
548 }
549
550 /**
551 * @brief Returns element wise Complex conjugate of Array
552 *
553 * @return Array<n, m, T>
554 */
555 inline Array<n, m, T> conj() const {
556 Array<n, m, T> res;
557 for (int i = 0; i < n * m; i++) {
558 res.c[i] = ::conj(c[i]);
559 }
560 return res;
561 }
562 /**
563 * @brief Returns real part of Array
564 *
565 * @return Array<n, m, T>
566 */
569 for (int i = 0; i < m * n; i++) {
570 res.c[i] = ::real(c[i]);
571 }
572 return res;
573 }
574
575 /// return imaginary part
578 for (int i = 0; i < m * n; i++) {
579 res.c[i] = ::imag(c[i]);
580 }
581 return res;
582 }
583
584 /// calculate square norm - sum of squared elements
585 hila::arithmetic_type<T> squarenorm() const {
586 hila::arithmetic_type<T> result = 0;
587 for (int i = 0; i < n * m; i++) {
588 result += ::squarenorm(c[i]);
589 }
590 return result;
591 }
592
593 /**
594 * @brief Fill Array with random elements
595 *
596 * @return Array<n, m, T>&
597 */
598 Array<n, m, T> &random() out_only {
599 for (int i = 0; i < n * m; i++) {
600 hila::random(c[i]);
601 }
602 return *this;
603 }
604
605 /**
606 * @brief Fill Array with Gaussian random elements
607 *
608 * @param width
609 * @return Array<n, m, T>&
610 */
611 Array<n, m, T> &gaussian_random(double width = 1.0) out_only {
612 for (int i = 0; i < n * m; i++) {
613 hila::gaussian_random(c[i], width);
614 }
615 return *this;
616 }
617
618 /// Convert to string for printing
619 std::string str(int prec = 8, char separator = ' ') const {
620 return this->asMatrix().str(prec, separator);
621 }
622
623 /// implement sort as casting to array
624#pragma hila novector
625 template <int N>
627 hila::sort order = hila::sort::ascending) const {
628 return this->asMatrix().sort(permutation, order).asArray();
629 }
630
631#pragma hila novector
632 Array<n, m, T> sort(hila::sort order = hila::sort::ascending) const {
633 return this->asMatrix().sort(order).asArray();
634 }
635};
636
637/**
638 * @brief Return conjugate Array
639 * @details Wrapper around Array::conj
640 *
641 * @tparam n Number of rows
642 * @tparam m Number of columns
643 * @tparam T Array element type
644 * @param arg Array to be conjugated
645 * @return Array<n, m, T>
646 */
647template <const int n, const int m, typename T>
649 return arg.conj();
650}
651/**
652 * @brief Return real part of Array
653 * @details Wrapper around Array::real
654 *
655 * @tparam n Number of rows
656 * @tparam m Number of columns
657 * @tparam T Array element type
658 * @param arg Array to return real part of
659 * @return Array<n, m, T>
660 */
661template <const int n, const int m, typename T>
663 return arg.real();
664}
665/**
666 * @brief Return imaginary part of Array
667 * @details Wrapper around Array::imag
668 *
669 * @tparam n Number of rows
670 * @tparam m Number of columns
671 * @tparam T Array element type
672 * @param arg Array to return real part of
673 * @return Array<n, m, T>
674 */
675template <const int n, const int m, typename T>
677 return arg.imag();
678}
679
680/**
681 * @brief Addition operator
682 * @details Defined for the following operations
683 *
684 * __Array + Array:__
685 *
686 * __NOTE__: Arrays must share same dimensionality
687 *
688 * \code {.cpp}
689 * Array<n,m,MyType> A, B, C;
690 * A = 1;
691 * B = 1;
692 * C = A + B; // C is uniformly 2
693 * \endcode
694 *
695 *
696 * __Scalar + Array / Array + Scalar:__
697 *
698 * __NOTE__: Exact definition exist in overloaded functions that can be viewed in source code.
699 * *
700 * \code {.cpp}
701 * Array<n,m,MyType> A,B;
702 * A = 1;
703 * B = A + 1; // B is uniformly 2
704 * \endcode
705 * @memberof Array
706 * @tparam n Number of rows
707 * @tparam m Number of columns
708 * @tparam T Array element type
709 * @param a
710 * @param b
711 * @return Array<n, m, T>
712 */
713template <int n, int m, typename T>
715 a += b;
716 return a;
717}
718
719/**
720 * @brief Subtraction operator
721 * @details Defined for the following operations
722 *
723 * __Array - Array:__
724 *
725 * __NOTE__: Arrays must share same dimensionality
726 *
727 * \code {.cpp}
728 * Array<n,m,MyType> A, B, C;
729 * A = 2;
730 * B = 1;
731 * C = A - B; // C is uniformly 2
732 * \endcode
733 *
734 *
735 * __Scalar - Array / Array - Scalar:__
736 *
737 * __NOTE__: Exact definition exist in overloaded functions that can be viewed in source code.
738 * *
739 * \code {.cpp}
740 * Array<n,m,MyType> A,B;
741 * A = 2;
742 * B = A - 1; // B is uniformly 2
743 * \endcode
744 * @memberof Array
745 * @tparam n Number of rows
746 * @tparam m Number of columns
747 * @tparam T Array element type
748 * @param a
749 * @param b
750 * @return Array<n, m, T>
751 */
752template <int n, int m, typename T>
754 a -= b;
755 return a;
756}
757
758// Array + scalar
759template <int n, int m, typename T, typename S,
760 std::enable_if_t<std::is_convertible<hila::type_plus<T, S>, T>::value, int> = 0>
761inline Array<n, m, T> operator+(Array<n, m, T> a, const S b) {
762 a += b;
763 return a;
764}
765
766// scalar + Array
767template <int n, int m, typename T, typename S,
768 std::enable_if_t<std::is_convertible<hila::type_plus<T, S>, T>::value, int> = 0>
769inline Array<n, m, T> operator+(const S b, Array<n, m, T> a) {
770 a += b;
771 return a;
772}
773
774// Array - scalar
775template <int n, int m, typename T, typename S,
776 std::enable_if_t<std::is_convertible<hila::type_minus<T, S>, T>::value, int> = 0>
777Array<n, m, T> operator-(Array<n, m, T> a, const S b) {
778 a -= b;
779 return a;
780}
781
782// scalar - Array
783template <int n, int m, typename T, typename S,
784 std::enable_if_t<std::is_convertible<hila::type_minus<S, T>, T>::value, int> = 0>
785inline Array<n, m, T> operator-(const S b, Array<n, m, T> a) {
786 for (int i = 0; i < n * m; i++)
787 a.c[i] = static_cast<T>(b) - a.c[i];
788 return a;
789}
790
791/**
792 * @brief Multiplication operator
793 * @details Defined for the following operations
794 *
795 * __Array*Array:__
796 *
797 * __NOTE__: Arrays must share same dimensionality
798 *
799 * \code {.cpp}
800 * Array<n,m,MyType> A, B, C;
801 * A = 2;
802 * B = 3;
803 * C = A*B; // C is uniformly 6
804 * \endcode
805 *
806 *
807 * __Scalar * Array / Array * Scalar:__
808 *
809 * __NOTE__: Exact definition exist in overloaded functions that can be viewed in source code.
810 * *
811 * \code {.cpp}
812 * Array<n,m,MyType> A,B;
813 * A = 2;
814 * B = A*3; // B is uniformly 6
815 * \endcode
816 * @memberof Array
817 * @tparam n Number of rows
818 * @tparam m Number of columns
819 * @tparam T Array element type
820 * @param a
821 * @param b
822 * @return Array<n, m, T>
823 */
824template <int n, int m, typename T>
826 a *= b;
827 return a;
828}
829
830/**
831 * @brief Division operator
832 * @details Defined for the following operations
833 *
834 * __Array/Array:__
835 *
836 * __NOTE__: Arrays must share same dimensionality
837 *
838 * \code {.cpp}
839 * Array<n,m,MyType> A, B, C;
840 * A = 4;
841 * B = 2;
842 * C = A/B; // C is uniformly 2
843 * \endcode
844 *
845 *
846 * __Scalar / Array / Array / Scalar:__
847 *
848 * __NOTE__: Exact definition exist in overloaded functions that can be viewed in source code.
849 * *
850 * \code {.cpp}
851 * Array<n,m,MyType> A,B;
852 * A = 4;
853 * B = A/2; // B is uniformly 2
854 * \endcode
855 * @memberof Array
856 * @tparam n Number of rows
857 * @tparam m Number of columns
858 * @tparam T Array element type
859 * @param a
860 * @param b
861 * @return Array<n, m, T>
862 */
863template <int n, int m, typename T>
865 a /= b;
866 return a;
867}
868
869// Array * scalar
870template <int n, int m, typename T, typename S,
871 std::enable_if_t<std::is_convertible<hila::type_mul<T, S>, T>::value, int> = 0>
872inline Array<n, m, T> operator*(Array<n, m, T> a, const S b) {
873 a *= b;
874 return a;
875}
876
877// scalar * Array
878template <int n, int m, typename T, typename S,
879 std::enable_if_t<std::is_convertible<hila::type_mul<T, S>, T>::value, int> = 0>
880inline Array<n, m, T> operator*(const S b, Array<n, m, T> a) {
881 a *= b;
882 return a;
883}
884
885// Array / scalar
886template <int n, int m, typename T, typename S,
887 std::enable_if_t<std::is_convertible<hila::type_div<T, S>, T>::value, int> = 0>
888inline Array<n, m, T> operator/(Array<n, m, T> a, const S b) {
889 a /= b;
890 return a;
891}
892
893// scalar / Array
894template <int n, int m, typename T, typename S,
895 std::enable_if_t<std::is_convertible<hila::type_div<S, T>, T>::value, int> = 0>
896inline Array<n, m, T> operator/(const S b, Array<n, m, T> a) {
897 for (int i = 0; i < n * m; i++)
898 a.c[i] = b / a.c[i];
899 return a;
900}
901
902
903/**
904 * @brief Stream operator
905 * @details Naive Stream operator for formatting Matrix for printing. Simply puts elements one after
906 * the other in row major order
907 *
908 * @tparam n Number of rows
909 * @tparam m Number of columns
910 * @tparam T Array element type
911 * @param strm OS stream
912 * @param A Array to write
913 * @return std::ostream&
914 */
915template <int n, int m, typename T>
916std::ostream &operator<<(std::ostream &strm, const Array<n, m, T> &A) {
917 return operator<<(strm, A.asMatrix());
918}
919
920namespace hila {
921
922/**
923 * @brief Converts Array object to string
924 *
925 * @tparam n Number of rows
926 * @tparam m Number of columns
927 * @tparam T Array element type
928 * @param A Array to convert to string
929 * @param prec Precision of T
930 * @param separator Separator between elements
931 * @return std::string
932 */
933template <int n, int m, typename T>
934std::string to_string(const Array<n, m, T> &A, int prec = 8, char separator = ' ') {
935 return to_string(A.asMatrix(), prec, separator);
936}
937
938template <int n, int m, typename T>
939std::string prettyprint(const Array<n, m, T> &A, int prec = 8) {
940 return prettyprint(A.asMatrix(), prec);
941}
942
943} // namespace hila
944
945
946/**
947 * @brief Return square norm of Array
948 * @details Wrapper around Array::squarenrom
949 *
950 * @tparam n Number of rows
951 * @tparam m Number of columns
952 * @tparam T Array element type
953 * @param rhs Array to compute squarenorm of
954 * @return hila::arithmetic_type<T>
955 */
956template <int n, int m, typename T>
957inline hila::arithmetic_type<T> squarenorm(const Array<n, m, T> &rhs) {
958 return rhs.squarenorm();
959}
960
961/** @name Arithmetic operations
962 * @details all operations are applied linearly to the Array A
963 *
964 * Most functions are self explanitory, but if necessary function will have a detailed section with
965 * additional information.
966 * @memberof Array
967 * @tparam n Number of rows
968 * @tparam m Number of columns
969 * @tparam T Array element type
970 * @param a Input array
971 * @return Array<n, m, T>
972 * @{
973 */
974
975/**
976 * @brief Square root
977 */
978template <int n, int m, typename T>
980 for (int i = 0; i < n * m; i++)
981 a.c[i] = sqrt(a.c[i]);
982 return a;
983}
984
985/**
986 * @brief Cuberoot
987 */
988template <int n, int m, typename T>
990 for (int i = 0; i < n * m; i++)
991 a.c[i] = cbrt(a.c[i]);
992 return a;
993}
994
995/**
996 * @brief Exponential
997 */
998template <int n, int m, typename T>
1000 for (int i = 0; i < n * m; i++)
1001 a.c[i] = exp(a.c[i]);
1002 return a;
1003}
1004
1005/**
1006 * @brief Logarithm
1007 */
1008template <int n, int m, typename T>
1010 for (int i = 0; i < n * m; i++)
1011 a.c[i] = log(a.c[i]);
1012 return a;
1013}
1014
1015/**
1016 * @brief Sine
1017 */
1018template <int n, int m, typename T>
1020 for (int i = 0; i < n * m; i++)
1021 a.c[i] = sin(a.c[i]);
1022 return a;
1023}
1024
1025/**
1026 * @brief Cosine
1027 */
1028template <int n, int m, typename T>
1030 for (int i = 0; i < n * m; i++)
1031 a.c[i] = cos(a.c[i]);
1032 return a;
1033}
1034
1035/**
1036 * @brief Tangent
1037 */
1038template <int n, int m, typename T>
1040 for (int i = 0; i < n * m; i++)
1041 a.c[i] = tan(a.c[i]);
1042 return a;
1043}
1044
1045/**
1046 * @brief Inverse Sine
1047 */
1048template <int n, int m, typename T>
1050 for (int i = 0; i < n * m; i++)
1051 a.c[i] = asin(a.c[i]);
1052 return a;
1053}
1054
1055/**
1056 * @brief Inverse Cosine
1057 */
1058template <int n, int m, typename T>
1060 for (int i = 0; i < n * m; i++)
1061 a.c[i] = acos(a.c[i]);
1062 return a;
1063}
1064
1065/**
1066 * @brief Inverse Tangent
1067 */
1068template <int n, int m, typename T>
1070 for (int i = 0; i < n * m; i++)
1071 a.c[i] = atan(a.c[i]);
1072 return a;
1073}
1074
1075/**
1076 * @brief Hyperbolic Sine
1077 */
1078template <int n, int m, typename T>
1080 for (int i = 0; i < n * m; i++)
1081 a.c[i] = sinh(a.c[i]);
1082 return a;
1083}
1084
1085/**
1086 * @brief Hyperbolic Cosine
1087 */
1088template <int n, int m, typename T>
1090 for (int i = 0; i < n * m; i++)
1091 a.c[i] = cosh(a.c[i]);
1092 return a;
1093}
1094
1095/**
1096 * @brief Hyperbolic tangent
1097 */
1098template <int n, int m, typename T>
1100 for (int i = 0; i < n * m; i++)
1101 a.c[i] = tanh(a.c[i]);
1102 return a;
1103}
1104
1105/**
1106 * @brief Inverse Hyperbolic Sine
1107 */
1108template <int n, int m, typename T>
1110 for (int i = 0; i < n * m; i++)
1111 a.c[i] = asinh(a.c[i]);
1112 return a;
1113}
1114
1115/**
1116 * @brief Inverse Hyperbolic Cosine
1117 */
1118template <int n, int m, typename T>
1120 for (int i = 0; i < n * m; i++)
1121 a.c[i] = acosh(a.c[i]);
1122 return a;
1123}
1124
1125/**
1126 * @brief Inverse Hyperbolic Tangent
1127 */
1128template <int n, int m, typename T>
1130 for (int i = 0; i < n * m; i++)
1131 a.c[i] = atanh(a.c[i]);
1132 return a;
1133}
1134
1135/**
1136 * @brief Power
1137 *
1138 * @details Array can be raised homogeneously to the power of a integer or scalar b or Array can be
1139 * raised to the power of another Array b element wise.
1140 *
1141 * @param b Array, Integer or Real scalar to raise to the power of
1142 */
1143template <int n, int m, typename T>
1145 for (int i = 0; i < n * m; i++)
1146 a.c[i] = pow(a.c[i], b);
1147 return a;
1148}
1149
1150template <int n, int m, typename T>
1151inline Array<n, m, T> pow(Array<n, m, T> a, T b) {
1152 for (int i = 0; i < n * m; i++)
1153 a.c[i] = pow(a.c[i], b);
1154 return a;
1155}
1156
1157template <int n, int m, typename T>
1158inline Array<n, m, T> pow(Array<n, m, T> a, const Array<n, m, T> &b) {
1159 for (int i = 0; i < n * m; i++)
1160 a.c[i] = pow(a.c[i], b.c[i]);
1161 return a;
1162}
1163
1164/**
1165 * @brief Rounding
1166 * @details Rounding function if T is arithmetic type
1167 * @tparam T Array element type, must be arithmetic type
1168 */
1169template <int n, int m, typename T, std::enable_if_t<hila::is_arithmetic<T>::value, int> = 0>
1171 for (int i = 0; i < n * m; i++)
1172 a.c[i] = round(a.c[i]);
1173 return a;
1174}
1175
1176/**
1177 * @brief Floor
1178 * @details Flooring function if T is arithmetic type
1179 * @tparam T Array element type, must be arithmetic type
1180 */
1181template <int n, int m, typename T, std::enable_if_t<hila::is_arithmetic<T>::value, int> = 0>
1183 for (int i = 0; i < n * m; i++)
1184 a.c[i] = floor(a.c[i]);
1185 return a;
1186}
1187
1188/**
1189 * @brief Ceiling
1190 * @details Ceiling function if T is arithmetic type
1191 * @tparam T Array element type, must be arithmetic type
1192 */
1193template <int n, int m, typename T, std::enable_if_t<hila::is_arithmetic<T>::value, int> = 0>
1195 for (int i = 0; i < n * m; i++)
1196 a.c[i] = ceil(a.c[i]);
1197 return a;
1198}
1199
1200/**
1201 * @brief Truncation
1202 * @details Truncation function if T is arithmetic type
1203 * @tparam T Array element type, must be arithmetic type
1204 */
1205template <int n, int m, typename T, std::enable_if_t<hila::is_arithmetic<T>::value, int> = 0>
1207 for (int i = 0; i < n * m; i++)
1208 a.c[i] = trunc(a.c[i]);
1209 return a;
1210}
1211/** @} */
1212
1213/**
1214 * @brief Array casting operation
1215 * @details Cast array to different number type or Complex type
1216 *
1217 * Allowed casting: number->number, number->Complex, Complex->Complex
1218 *
1219 * Not allowed casting: Complex->number
1220 * @tparam Ntype Type to cast to
1221 * @tparam T Input Array type
1222 * @tparam n Number of rows
1223 * @tparam m Number of columns
1224 * @param mat Array to cast
1225 * @return Array<n, m, Ntype>
1226 */
1227template <typename Ntype, typename T, int n, int m,
1228 std::enable_if_t<hila::is_arithmetic<T>::value, int> = 0>
1231 for (int i = 0; i < n * m; i++)
1232 res.c[i] = mat.c[i];
1233 return res;
1234}
1235
1236template <typename Ntype, typename T, int n, int m,
1237 std::enable_if_t<hila::is_complex<T>::value, int> = 0>
1238Array<n, m, Ntype> cast_to(const Array<n, m, T> &mat) {
1240 for (int i = 0; i < n * m; i++)
1241 res.c[i] = cast_to<Ntype>(mat.c[i]);
1242 return res;
1243}
1244
1245
1246/// @brief Array1d is alias to Array where m = 1
1247template <int n, typename T = double>
1249
1250/// @brief Array2d is simply alias of Array
1251template <int n, int m, typename T = double>
1253
1254#endif
Array< n, m, T > conj(const Array< n, m, T > &arg)
Return conjugate Array.
Definition array.h:648
std::ostream & operator<<(std::ostream &strm, const Array< n, m, T > &A)
Stream operator.
Definition array.h:916
Array< n, m, hila::arithmetic_type< T > > imag(const Array< n, m, T > &arg)
Return imaginary part of Array.
Definition array.h:676
hila::arithmetic_type< T > squarenorm(const Array< n, m, T > &rhs)
Return square norm of Array.
Definition array.h:957
Array< n, m, hila::arithmetic_type< T > > real(const Array< n, m, T > &arg)
Return real part of Array.
Definition array.h:662
Array type
Definition array.h:43
Array< n, m, T > operator/(Array< n, m, T > a, const Array< n, m, T > &b)
Division operator.
Definition array.h:864
Array(const S rhs)
Scalar constructor.
Definition array.h:100
Array(std::initializer_list< S > rhs)
Initializer list constructor.
Definition array.h:134
Array< n, m, T > & operator+=(const Array< n, m, S > &rhs)
Add assign operator with Array or scalar.
Definition array.h:400
Array< n, m, T > floor(Array< n, m, T > a)
Floor.
Definition array.h:1182
constexpr int rows() const
Returns number of rows.
Definition array.h:147
Array(const std::nullptr_t &z)
Zero constructor.
Definition array.h:113
constexpr int size() const
Returns size of Vector Array or square Array.
Definition array.h:167
T & e(const int i, const int j)
Const overload.
Definition array.h:203
std::string str(int prec=8, char separator=' ') const
Convert to string for printing.
Definition array.h:619
Array< n, m, T > round(Array< n, m, T > a)
Rounding.
Definition array.h:1170
Array< n, m, T > operator+() const
Unary + operator.
Definition array.h:300
Array< n, m, T > asin(Array< n, m, T > a)
Inverse Sine.
Definition array.h:1049
Array< n, m, T > & operator*=(const Array< n, m, S > &rhs)
Multiply assign scalar or array.
Definition array.h:488
bool operator==(const Array< n1, m1, S > &rhs) const
Compare equality of Arrays.
Definition array.h:340
Array< n, m, T > exp(Array< n, m, T > a)
Exponential.
Definition array.h:999
T e(const int i) const
Const overload.
Definition array.h:210
Array< n, m, T > tanh(Array< n, m, T > a)
Hyperbolic tangent.
Definition array.h:1099
Array< n, m, T > trunc(Array< n, m, T > a)
Truncation.
Definition array.h:1206
T & e(const int i)
Standard array indexing operation for 1D Array.
Definition array.h:225
Array< n, m, T > sort(Vector< N, int > &permutation, hila::sort order=hila::sort::ascending) const
implement sort as casting to array
Definition array.h:626
Array< n, m, hila::arithmetic_type< T > > real() const
Returns real part of Array.
Definition array.h:567
Array< n, m, T > cosh(Array< n, m, T > a)
Hyperbolic Cosine.
Definition array.h:1089
Array< n, m, T > operator-(Array< n, m, T > a, const Array< n, m, T > &b)
Subtraction operator.
Definition array.h:753
Array< n, m, T > tan(Array< n, m, T > a)
Tangent.
Definition array.h:1039
Array()=default
Default constructor.
Array< n, m, hila::arithmetic_type< T > > imag() const
return imaginary part
Definition array.h:576
Array< n, m, T > & operator/=(const Array< n, m, S > &rhs)
Division assign with array or scalar.
Definition array.h:533
Array< n, m, T > operator-() const
Unary - operator.
Definition array.h:286
constexpr int columns() const
Returns number of columns.
Definition array.h:155
Array< n, m, T > operator+(Array< n, m, T > a, const Array< n, m, T > &b)
Addition operator.
Definition array.h:714
T e(const int i, const int j) const
Standard array indexing operation for 2D Array.
Definition array.h:197
Array< n, m, T > cos(Array< n, m, T > a)
Cosine.
Definition array.h:1029
bool operator!=(const Array< n1, m1, S > &rhs) const
Compare non-equality of two Arrays.
Definition array.h:366
Array< n, m, T > cbrt(Array< n, m, T > a)
Cuberoot.
Definition array.h:989
hila::arithmetic_type< T > squarenorm() const
calculate square norm - sum of squared elements
Definition array.h:585
Vector< n, T > & asVector()
Cast Array1D to Vector.
Definition array.h:248
Array< n, m, T > acos(Array< n, m, T > a)
Inverse Cosine.
Definition array.h:1059
Array< n, m, T > & operator*=(const S rhs)
multiply assign with scalar
Definition array.h:498
Array< n, m, T > acosh(Array< n, m, T > a)
Inverse Hyperbolic Cosine.
Definition array.h:1119
Array< n, m, T > sqrt(Array< n, m, T > a)
Square root.
Definition array.h:979
Array< n, m, T > conj() const
Returns element wise Complex conjugate of Array.
Definition array.h:555
Array< n, m, T > asinh(Array< n, m, T > a)
Inverse Hyperbolic Sine.
Definition array.h:1109
Matrix< n, m, T > & asMatrix()
Cast Array to Matrix.
Definition array.h:234
Array< n, m, T > & random()
Fill Array with random elements.
Definition array.h:598
Array< n, m, T > ceil(Array< n, m, T > a)
Ceiling.
Definition array.h:1194
Array< n, m, T > atan(Array< n, m, T > a)
Inverse Tangent.
Definition array.h:1069
Array< n, m, T > pow(Array< n, m, T > a, int b)
Power.
Definition array.h:1144
Array< n, m, T > & operator-=(const Array< n, m, S > &rhs)
Subtract assign operator with Array or scalar.
Definition array.h:433
Array(const Array< n, m, T > &v)=default
Copy constructor.
Array< n, m, T > operator*(Array< n, m, T > a, const Array< n, m, T > &b)
Multiplication operator.
Definition array.h:825
Array< n, m, Ntype > cast_to(const Array< n, m, T > &mat)
Array casting operation.
Definition array.h:1229
Array< n, m, T > sin(Array< n, m, T > a)
Sine.
Definition array.h:1019
Array< n, m, T > sinh(Array< n, m, T > a)
Hyperbolic Sine.
Definition array.h:1079
Array< n, m, T > & gaussian_random(double width=1.0)
Fill Array with Gaussian random elements.
Definition array.h:611
Array< n, m, T > atanh(Array< n, m, T > a)
Inverse Hyperbolic Tangent.
Definition array.h:1129
Array< n, m, T > & operator=(const S rhs)
Scalar assignment operator.
Definition array.h:321
Array< n, m, T > log(Array< n, m, T > a)
Logarithm.
Definition array.h:1009
Define type DiagonalMatrix<n,T>
Matrix class which defines matrix operations.
Definition matrix.h:1679
T arg(const Complex< T > &a)
Return argument of Complex number.
Definition cmplx.h:1334
Definition of Matrix types.
Invert diagonal + const. matrix using Sherman-Morrison formula.
Definition array.h:920
T gaussian_random()
Template function T hila::gaussian_random<T>(),generates gaussian random value of type T,...
Definition random.h:183
double random()
Real valued uniform random number generator.
Definition hila_gpu.cpp:120
std::string to_string(const Array< n, m, T > &A, int prec=8, char separator=' ')
Converts Array object to string.
Definition array.h:934
hila::is_complex_or_arithmetic<T>::value
Definition cmplx.h:750