Object

org.ekrich.blas.snic

blas

Related Doc: package snic

Permalink

object blas

Scala Native extern C interface to CBLAS Version 3.8.0

Annotations
@link( "cblas" ) @extern()
Linear Supertypes
AnyRef, Any
Ordering
  1. Alphabetic
  2. By Inheritance
Inherited
  1. blas
  2. AnyRef
  3. Any
  1. Hide All
  2. Show All
Visibility
  1. Public
  2. All

Type Members

  1. type CblasIndex = Long

    Permalink

Value Members

  1. final def !=(arg0: Any): Boolean

    Permalink
    Definition Classes
    AnyRef → Any
  2. final def ##(): Int

    Permalink
    Definition Classes
    AnyRef → Any
  3. final def ==(arg0: Any): Boolean

    Permalink
    Definition Classes
    AnyRef → Any
  4. final def asInstanceOf[T0]: T0

    Permalink
    Definition Classes
    Any
  5. def cblas_caxpy(N: CInt, alpha: Ptr[CFloatComplex], X: Ptr[CFloatComplex], incX: CInt, Y: Ptr[CFloatComplex], incY: CInt): Unit

    Permalink

    A constant times a vector plus a vector (single-precision complex).

  6. def cblas_ccopy(N: CInt, X: Ptr[CFloatComplex], incX: CInt, Y: Ptr[CFloatComplex], incY: CInt): Unit

    Permalink

    Copies a vector to another vector (single-precision complex).

  7. def cblas_cdotc_sub(N: CInt, X: Ptr[CFloatComplex], incX: CInt, Y: Ptr[CFloatComplex], incY: CInt, dotc: Ptr[CFloatComplex]): Unit

    Permalink

    Dot product of the complex conjugate of a single-precision complex vector with a second single-precision complex vector.

    Dot product of the complex conjugate of a single-precision complex vector with a second single-precision complex vector.

    dotc

    - The result vector. Computes conjg(X) * Y.

  8. def cblas_cdotu_sub(N: CInt, X: Ptr[CFloatComplex], incX: CInt, Y: Ptr[CFloatComplex], incY: CInt, dotu: Ptr[CFloatComplex]): Unit

    Permalink

    Dot product of two single-precision complex vectors.

    Dot product of two single-precision complex vectors.

    dotu

    - The result vector.

  9. def cblas_cgbmv(order: CBLAS_ORDER, TransA: CBLAS_TRANSPOSE, M: CInt, N: CInt, KL: CInt, KU: CInt, alpha: Ptr[CFloatComplex], A: Ptr[CFloatComplex], lda: CInt, X: Ptr[CFloatComplex], incX: CInt, beta: Ptr[CFloatComplex], Y: Ptr[CFloatComplex], incY: CInt): Unit

    Permalink
  10. def cblas_cgemm(Order: CBLAS_ORDER, TransA: CBLAS_TRANSPOSE, TransB: CBLAS_TRANSPOSE, M: CInt, N: CInt, K: CInt, alpha: Ptr[CFloatComplex], A: Ptr[CFloatComplex], lda: CInt, B: Ptr[CFloatComplex], ldb: CInt, beta: Ptr[CFloatComplex], C: Ptr[CFloatComplex], ldc: CInt): Unit

    Permalink
  11. def cblas_cgemv(order: CBLAS_ORDER, TransA: CBLAS_TRANSPOSE, M: CInt, N: CInt, alpha: Ptr[CFloatComplex], A: Ptr[CFloatComplex], lda: CInt, X: Ptr[CFloatComplex], incX: CInt, beta: Ptr[CFloatComplex], Y: Ptr[CFloatComplex], incY: CInt): Unit

    Permalink
  12. def cblas_cgerc(order: CBLAS_ORDER, M: CInt, N: CInt, alpha: Ptr[CFloatComplex], X: Ptr[CFloatComplex], incX: CInt, Y: Ptr[CFloatComplex], incY: CInt, A: Ptr[CFloatComplex], lda: CInt): Unit

    Permalink
  13. def cblas_cgeru(order: CBLAS_ORDER, M: CInt, N: CInt, alpha: Ptr[CFloatComplex], X: Ptr[CFloatComplex], incX: CInt, Y: Ptr[CFloatComplex], incY: CInt, A: Ptr[CFloatComplex], lda: CInt): Unit

    Permalink
  14. def cblas_chbmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, K: CInt, alpha: Ptr[CFloatComplex], A: Ptr[CFloatComplex], lda: CInt, X: Ptr[CFloatComplex], incX: CInt, beta: Ptr[CFloatComplex], Y: Ptr[CFloatComplex], incY: CInt): Unit

    Permalink
  15. def cblas_chemm(Order: CBLAS_ORDER, Side: CBLAS_SIDE, Uplo: CBLAS_UPLO, M: CInt, N: CInt, alpha: Ptr[CFloatComplex], A: Ptr[CFloatComplex], lda: CInt, B: Ptr[CFloatComplex], ldb: CInt, beta: Ptr[CFloatComplex], C: Ptr[CFloatComplex], ldc: CInt): Unit

    Permalink
  16. def cblas_chemv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: Ptr[CFloatComplex], A: Ptr[CFloatComplex], lda: CInt, X: Ptr[CFloatComplex], incX: CInt, beta: Ptr[CFloatComplex], Y: Ptr[CFloatComplex], incY: CInt): Unit

    Permalink
  17. def cblas_cher(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: CFloat, X: Ptr[CFloatComplex], incX: CInt, A: Ptr[CFloatComplex], lda: CInt): Unit

    Permalink
  18. def cblas_cher2(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: Ptr[CFloatComplex], X: Ptr[CFloatComplex], incX: CInt, Y: Ptr[CFloatComplex], incY: CInt, A: Ptr[CFloatComplex], lda: CInt): Unit

    Permalink
  19. def cblas_cher2k(Order: CBLAS_ORDER, Uplo: CBLAS_UPLO, Trans: CBLAS_TRANSPOSE, N: CInt, K: CInt, alpha: Ptr[CFloatComplex], A: Ptr[CFloatComplex], lda: CInt, B: Ptr[CFloatComplex], ldb: CInt, beta: CFloat, C: Ptr[CFloatComplex], ldc: CInt): Unit

    Permalink
  20. def cblas_cherk(Order: CBLAS_ORDER, Uplo: CBLAS_UPLO, Trans: CBLAS_TRANSPOSE, N: CInt, K: CInt, alpha: CFloat, A: Ptr[CFloatComplex], lda: CInt, beta: CFloat, C: Ptr[CFloatComplex], ldc: CInt): Unit

    Permalink
  21. def cblas_chpmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: Ptr[CFloatComplex], Ap: Ptr[CFloatComplex], X: Ptr[CFloatComplex], incX: CInt, beta: Ptr[CFloatComplex], Y: Ptr[CFloatComplex], incY: CInt): Unit

    Permalink
  22. def cblas_chpr(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: CFloat, X: Ptr[CFloatComplex], incX: CInt, A: Ptr[CFloatComplex]): Unit

    Permalink
  23. def cblas_chpr2(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: Ptr[CFloatComplex], X: Ptr[CFloatComplex], incX: CInt, Y: Ptr[CFloatComplex], incY: CInt, Ap: Ptr[CFloatComplex]): Unit

    Permalink
  24. def cblas_cscal(N: CInt, alpha: Ptr[CFloatComplex], X: Ptr[CFloatComplex], incX: CInt): Unit

    Permalink

    Multiplies each element of a vector by a constant (single-precision complex).

  25. def cblas_csscal(N: CInt, alpha: CFloat, X: Ptr[CFloatComplex], incX: CInt): Unit

    Permalink

    Multiplies each element of a vector by a constant (single-precision complex).

  26. def cblas_cswap(N: CInt, X: Ptr[CFloatComplex], incX: CInt, Y: Ptr[CFloatComplex], incY: CInt): Unit

    Permalink

    Exchanges the elements of two vectors (single-precision complex).

  27. def cblas_csymm(Order: CBLAS_ORDER, Side: CBLAS_SIDE, Uplo: CBLAS_UPLO, M: CInt, N: CInt, alpha: Ptr[CFloatComplex], A: Ptr[CFloatComplex], lda: CInt, B: Ptr[CFloatComplex], ldb: CInt, beta: Ptr[CFloatComplex], C: Ptr[CFloatComplex], ldc: CInt): Unit

    Permalink
  28. def cblas_csyr2k(Order: CBLAS_ORDER, Uplo: CBLAS_UPLO, Trans: CBLAS_TRANSPOSE, N: CInt, K: CInt, alpha: Ptr[CFloatComplex], A: Ptr[CFloatComplex], lda: CInt, B: Ptr[CFloatComplex], ldb: CInt, beta: Ptr[CFloatComplex], C: Ptr[CFloatComplex], ldc: CInt): Unit

    Permalink
  29. def cblas_csyrk(Order: CBLAS_ORDER, Uplo: CBLAS_UPLO, Trans: CBLAS_TRANSPOSE, N: CInt, K: CInt, alpha: Ptr[CFloatComplex], A: Ptr[CFloatComplex], lda: CInt, beta: Ptr[CFloatComplex], C: Ptr[CFloatComplex], ldc: CInt): Unit

    Permalink
  30. def cblas_ctbmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, K: CInt, A: Ptr[CFloatComplex], lda: CInt, X: Ptr[CFloatComplex], incX: CInt): Unit

    Permalink
  31. def cblas_ctbsv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, K: CInt, A: Ptr[CFloatComplex], lda: CInt, X: Ptr[CFloatComplex], incX: CInt): Unit

    Permalink
  32. def cblas_ctpmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, Ap: Ptr[CFloatComplex], X: Ptr[CFloatComplex], incX: CInt): Unit

    Permalink
  33. def cblas_ctpsv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, Ap: Ptr[CFloatComplex], X: Ptr[CFloatComplex], incX: CInt): Unit

    Permalink
  34. def cblas_ctrmm(Order: CBLAS_ORDER, Side: CBLAS_SIDE, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, M: CInt, N: CInt, alpha: Ptr[CFloatComplex], A: Ptr[CFloatComplex], lda: CInt, B: Ptr[CFloatComplex], ldb: CInt): Unit

    Permalink
  35. def cblas_ctrmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, A: Ptr[CFloatComplex], lda: CInt, X: Ptr[CFloatComplex], incX: CInt): Unit

    Permalink
  36. def cblas_ctrsm(Order: CBLAS_ORDER, Side: CBLAS_SIDE, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, M: CInt, N: CInt, alpha: Ptr[CFloatComplex], A: Ptr[CFloatComplex], lda: CInt, B: Ptr[CFloatComplex], ldb: CInt): Unit

    Permalink
  37. def cblas_ctrsv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, A: Ptr[CFloatComplex], lda: CInt, X: Ptr[CFloatComplex], incX: CInt): Unit

    Permalink
  38. def cblas_dasum(N: CInt, X: Ptr[CDouble], incX: CInt): CDouble

    Permalink

    Sum of the absolute values of elements in a vector (double-precision).

  39. def cblas_daxpy(N: CInt, alpha: CDouble, X: Ptr[CDouble], incX: CInt, Y: Ptr[CDouble], incY: CInt): Unit

    Permalink

    Computes a constant times a vector plus a vector (double-precision).

  40. def cblas_dcopy(N: CInt, X: Ptr[CDouble], incX: CInt, Y: Ptr[CDouble], incY: CInt): Unit

    Permalink

    Copies a vector to another vector (double-precision).

  41. def cblas_ddot(N: CInt, X: Ptr[CDouble], incX: CInt, Y: Ptr[CDouble], incY: CInt): CDouble

    Permalink

    Dot product of two vectors (double-precision).

  42. def cblas_dgbmv(order: CBLAS_ORDER, TransA: CBLAS_TRANSPOSE, M: CInt, N: CInt, KL: CInt, KU: CInt, alpha: CDouble, A: Ptr[CDouble], lda: CInt, X: Ptr[CDouble], incX: CInt, beta: CDouble, Y: Ptr[CDouble], incY: CInt): Unit

    Permalink
  43. def cblas_dgemm(Order: CBLAS_ORDER, TransA: CBLAS_TRANSPOSE, TransB: CBLAS_TRANSPOSE, M: CInt, N: CInt, K: CInt, alpha: CDouble, A: Ptr[CDouble], lda: CInt, B: Ptr[CDouble], ldb: CInt, beta: CDouble, C: Ptr[CDouble], ldc: CInt): Unit

    Permalink
  44. def cblas_dgemv(order: CBLAS_ORDER, TransA: CBLAS_TRANSPOSE, M: CInt, N: CInt, alpha: CDouble, A: Ptr[CDouble], lda: CInt, X: Ptr[CDouble], incX: CInt, beta: CDouble, Y: Ptr[CDouble], incY: CInt): Unit

    Permalink
  45. def cblas_dger(order: CBLAS_ORDER, M: CInt, N: CInt, alpha: CDouble, X: Ptr[CDouble], incX: CInt, Y: Ptr[CDouble], incY: CInt, A: Ptr[CDouble], lda: CInt): Unit

    Permalink
  46. def cblas_dnrm2(N: CInt, X: Ptr[CDouble], incX: CInt): CDouble

    Permalink

    L2 norm (Euclidian length) of a vector (double-precision).

  47. def cblas_drot(N: CInt, X: Ptr[CDouble], incX: CInt, Y: Ptr[CDouble], incY: CInt, c: CDouble, s: CDouble): Unit

    Permalink

    Applies a Givens rotation matrix to a pair of vectors.

  48. def cblas_drotg(a: Ptr[CDouble], b: Ptr[CDouble], c: Ptr[CDouble], s: Ptr[CDouble]): Unit

    Permalink

    Constructs a Givens rotation matrix.

  49. def cblas_drotm(N: CInt, X: Ptr[CDouble], incX: CInt, Y: Ptr[CDouble], incY: CInt, P: Ptr[CDouble]): Unit

    Permalink

    Applies a modified Givens transformation (single precision).

  50. def cblas_drotmg(d1: Ptr[CDouble], d2: Ptr[CDouble], b1: Ptr[CDouble], b2: CDouble, P: Ptr[CDouble]): Unit

    Permalink

    Generates a modified Givens rotation matrix.

  51. def cblas_dsbmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, K: CInt, alpha: CDouble, A: Ptr[CDouble], lda: CInt, X: Ptr[CDouble], incX: CInt, beta: CDouble, Y: Ptr[CDouble], incY: CInt): Unit

    Permalink
  52. def cblas_dscal(N: CInt, alpha: CDouble, X: Ptr[CDouble], incX: CInt): Unit

    Permalink

    Multiplies each element of a vector by a constant (double-precision).

  53. def cblas_dsdot(N: CInt, X: Ptr[CFloat], incX: CInt, Y: Ptr[CFloat], incY: CInt): CDouble

    Permalink

    Double-precision dot product of a pair of single-precision vectors.

  54. def cblas_dspmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: CDouble, Ap: Ptr[CDouble], X: Ptr[CDouble], incX: CInt, beta: CDouble, Y: Ptr[CDouble], incY: CInt): Unit

    Permalink
  55. def cblas_dspr(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: CDouble, X: Ptr[CDouble], incX: CInt, Ap: Ptr[CDouble]): Unit

    Permalink
  56. def cblas_dspr2(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: CDouble, X: Ptr[CDouble], incX: CInt, Y: Ptr[CDouble], incY: CInt, A: Ptr[CDouble]): Unit

    Permalink
  57. def cblas_dswap(N: CInt, X: Ptr[CDouble], incX: CInt, Y: Ptr[CDouble], incY: CInt): Unit

    Permalink

    Exchanges the elements of two vectors (double precision).

  58. def cblas_dsymm(Order: CBLAS_ORDER, Side: CBLAS_SIDE, Uplo: CBLAS_UPLO, M: CInt, N: CInt, alpha: CDouble, A: Ptr[CDouble], lda: CInt, B: Ptr[CDouble], ldb: CInt, beta: CDouble, C: Ptr[CDouble], ldc: CInt): Unit

    Permalink
  59. def cblas_dsymv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: CDouble, A: Ptr[CDouble], lda: CInt, X: Ptr[CDouble], incX: CInt, beta: CDouble, Y: Ptr[CDouble], incY: CInt): Unit

    Permalink
  60. def cblas_dsyr(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: CDouble, X: Ptr[CDouble], incX: CInt, A: Ptr[CDouble], lda: CInt): Unit

    Permalink
  61. def cblas_dsyr2(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: CDouble, X: Ptr[CDouble], incX: CInt, Y: Ptr[CDouble], incY: CInt, A: Ptr[CDouble], lda: CInt): Unit

    Permalink
  62. def cblas_dsyr2k(Order: CBLAS_ORDER, Uplo: CBLAS_UPLO, Trans: CBLAS_TRANSPOSE, N: CInt, K: CInt, alpha: CDouble, A: Ptr[CDouble], lda: CInt, B: Ptr[CDouble], ldb: CInt, beta: CDouble, C: Ptr[CDouble], ldc: CInt): Unit

    Permalink
  63. def cblas_dsyrk(Order: CBLAS_ORDER, Uplo: CBLAS_UPLO, Trans: CBLAS_TRANSPOSE, N: CInt, K: CInt, alpha: CDouble, A: Ptr[CDouble], lda: CInt, beta: CDouble, C: Ptr[CDouble], ldc: CInt): Unit

    Permalink
  64. def cblas_dtbmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, K: CInt, A: Ptr[CDouble], lda: CInt, X: Ptr[CDouble], incX: CInt): Unit

    Permalink
  65. def cblas_dtbsv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, K: CInt, A: Ptr[CDouble], lda: CInt, X: Ptr[CDouble], incX: CInt): Unit

    Permalink
  66. def cblas_dtpmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, Ap: Ptr[CDouble], X: Ptr[CDouble], incX: CInt): Unit

    Permalink
  67. def cblas_dtpsv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, Ap: Ptr[CDouble], X: Ptr[CDouble], incX: CInt): Unit

    Permalink
  68. def cblas_dtrmm(Order: CBLAS_ORDER, Side: CBLAS_SIDE, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, M: CInt, N: CInt, alpha: CDouble, A: Ptr[CDouble], lda: CInt, B: Ptr[CDouble], ldb: CInt): Unit

    Permalink
  69. def cblas_dtrmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, A: Ptr[CDouble], lda: CInt, X: Ptr[CDouble], incX: CInt): Unit

    Permalink
  70. def cblas_dtrsm(Order: CBLAS_ORDER, Side: CBLAS_SIDE, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, M: CInt, N: CInt, alpha: CDouble, A: Ptr[CDouble], lda: CInt, B: Ptr[CDouble], ldb: CInt): Unit

    Permalink
  71. def cblas_dtrsv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, A: Ptr[CDouble], lda: CInt, X: Ptr[CDouble], incX: CInt): Unit

    Permalink
  72. def cblas_dzasum(N: CInt, X: Ptr[CDoubleComplex], incX: CInt): CDouble

    Permalink

    Sum of the absolute values of real and imaginary parts of elements in a vector (double-precision complex).

  73. def cblas_dznrm2(N: CInt, X: Ptr[CDoubleComplex], incX: CInt): CDouble

    Permalink

    Unitary norm of a vector (double-precision complex).

  74. def cblas_icamax(N: CInt, X: Ptr[CFloatComplex], incX: CInt): CblasIndex

    Permalink

    returns

    Index of the element with the largest absolute value in a vector (single-precision complex).

  75. def cblas_idamax(N: CInt, X: Ptr[CDouble], incX: CInt): CblasIndex

    Permalink

    returns

    Index of the element with the largest absolute value in a vector (double-precision).

  76. def cblas_isamax(N: CInt, X: Ptr[CFloat], incX: CInt): CblasIndex

    Permalink

    returns

    Index of the element with the largest absolute value in a vector (single-precision).

  77. def cblas_izamax(N: CInt, X: Ptr[CDoubleComplex], incX: CInt): CblasIndex

    Permalink

    returns

    Index of the element with the largest absolute value in a vector (double-precision complex).

  78. def cblas_sasum(N: CInt, X: Ptr[CFloat], incX: CInt): CFloat

    Permalink

    Sum of the absolute values of elements in a vector (single-precision).

  79. def cblas_saxpy(N: CInt, alpha: CFloat, X: Ptr[CFloat], incX: CInt, Y: Ptr[CFloat], incY: CInt): Unit

    Permalink

    A constant times a vector plus a vector (single-precision).

    A constant times a vector plus a vector (single-precision).

    alpha

    The initial value to add to the dot product.

  80. def cblas_scasum(N: CInt, X: Ptr[CFloatComplex], incX: CInt): CFloat

    Permalink

    Sum of the absolute values of real and imaginary parts of elements in a vector (single-precision complex).

  81. def cblas_scnrm2(N: CInt, X: Ptr[CFloatComplex], incX: CInt): CFloat

    Permalink

    Unitary norm of a vector (single-precision complex).

  82. def cblas_scopy(N: CInt, X: Ptr[CFloat], incX: CInt, Y: Ptr[CFloat], incY: CInt): Unit

    Permalink

    Copies a vector to another vector (single-precision).

  83. def cblas_sdot(N: CInt, X: Ptr[CFloat], incX: CInt, Y: Ptr[CFloat], incY: CInt): CFloat

    Permalink

    Dot product of two vectors (single-precision).

  84. def cblas_sdsdot(N: CInt, alpha: CFloat, X: Ptr[CFloat], incX: CInt, Y: Ptr[CFloat], incY: CInt): CFloat

    Permalink

    Dot product of two single-precision vectors plus an initial single-precision value.

    Dot product of two single-precision vectors plus an initial single-precision value.

    N

    The number of elements in the vectors.

    alpha

    The initial value to add to the dot product.

    X

    Vector X.

    incX

    Stride within X. For example, if incX is 7, every 7th element is used.

    Y

    Vector Y.

    incY

    Stride within Y. For example, if incY is 7, every 7th element is used.

    returns

    See description above.

  85. def cblas_sgbmv(order: CBLAS_ORDER, TransA: CBLAS_TRANSPOSE, M: CInt, N: CInt, KL: CInt, KU: CInt, alpha: CFloat, A: Ptr[CFloat], lda: CInt, X: Ptr[CFloat], incX: CInt, beta: CFloat, Y: Ptr[CFloat], incY: CInt): Unit

    Permalink
  86. def cblas_sgemm(Order: CBLAS_ORDER, TransA: CBLAS_TRANSPOSE, TransB: CBLAS_TRANSPOSE, M: CInt, N: CInt, K: CInt, alpha: CFloat, A: Ptr[CFloat], lda: CInt, B: Ptr[CFloat], ldb: CInt, beta: CFloat, C: Ptr[CFloat], ldc: CInt): Unit

    Permalink
  87. def cblas_sgemv(order: CBLAS_ORDER, TransA: CBLAS_TRANSPOSE, M: CInt, N: CInt, alpha: CFloat, A: Ptr[CFloat], lda: CInt, X: Ptr[CFloat], incX: CInt, beta: CFloat, Y: Ptr[CFloat], incY: CInt): Unit

    Permalink
  88. def cblas_sger(order: CBLAS_ORDER, M: CInt, N: CInt, alpha: CFloat, X: Ptr[CFloat], incX: CInt, Y: Ptr[CFloat], incY: CInt, A: Ptr[CFloat], lda: CInt): Unit

    Permalink
  89. def cblas_snrm2(N: CInt, X: Ptr[CFloat], incX: CInt): CFloat

    Permalink

    L2 norm (Euclidian length) of a vector (single-precision).

  90. def cblas_srot(N: CInt, X: Ptr[CFloat], incX: CInt, Y: Ptr[CFloat], incY: CInt, c: CFloat, s: CFloat): Unit

    Permalink

    Applies a Givens rotation matrix to a pair of vectors.

  91. def cblas_srotg(a: Ptr[CFloat], b: Ptr[CFloat], c: Ptr[CFloat], s: Ptr[CFloat]): Unit

    Permalink

    Constructs a Givens rotation matrix.

    Constructs a Givens rotation matrix.

    a

    Single-precision value a. Overwritten on return with result r.

    b

    Single-precision value b. Overwritten on return with result z (zero).

    c

    Unused on entry. Overwritten on return with the value cos(θ).

    s

    Unused on entry. Overwritten on return with the value sin(θ).

  92. def cblas_srotm(N: CInt, X: Ptr[CFloat], incX: CInt, Y: Ptr[CFloat], incY: CInt, P: Ptr[CFloat]): Unit

    Permalink

    Applies a modified Givens transformation (single precision).

  93. def cblas_srotmg(d1: Ptr[CFloat], d2: Ptr[CFloat], b1: Ptr[CFloat], b2: CFloat, P: Ptr[CFloat]): Unit

    Permalink

    Generates a modified Givens rotation matrix.

    Generates a modified Givens rotation matrix.

    d1

    Scaling factor D1.

    d2

    Scaling factor D2.

    b1

    Scaling factor B1.

    b2

    Scaling factor B2.

    P

    A 5-element vector: P[0] Flag value that defines the form of matrix H. -2.0: matrix H contains the identity matrix. -1.0: matrix H is identical to matrix SH (defined by the remaining values in the vector). 0.0: H[1,2] and H[2,1] are obtained from matrix SH. The remaining values are both 1.0. 1.0: H[1,1] and H[2,2] are obtained from matrix SH. H[1,2] is 1.0. H[2,1] is -1.0. P[1] Contains SH[1,1]. P[2] Contains SH[2,1]. P[3] Contains SH[1,2]. P[4] Contains SH[2,2].

  94. def cblas_ssbmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, K: CInt, alpha: CFloat, A: Ptr[CFloat], lda: CInt, X: Ptr[CFloat], incX: CInt, beta: CFloat, Y: Ptr[CFloat], incY: CInt): Unit

    Permalink
  95. def cblas_sscal(N: CInt, alpha: CFloat, X: Ptr[CFloat], incX: CInt): Unit

    Permalink

    Multiplies each element of a vector by a constant (single-precision).

  96. def cblas_sspmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: CFloat, Ap: Ptr[CFloat], X: Ptr[CFloat], incX: CInt, beta: CFloat, Y: Ptr[CFloat], incY: CInt): Unit

    Permalink
  97. def cblas_sspr(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: CFloat, X: Ptr[CFloat], incX: CInt, Ap: Ptr[CFloat]): Unit

    Permalink
  98. def cblas_sspr2(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: CFloat, X: Ptr[CFloat], incX: CInt, Y: Ptr[CFloat], incY: CInt, A: Ptr[CFloat]): Unit

    Permalink
  99. def cblas_sswap(N: CInt, X: Ptr[CFloat], incX: CInt, Y: Ptr[CFloat], incY: CInt): Unit

    Permalink

    Exchanges the elements of two vectors (single precision).

    Exchanges the elements of two vectors (single precision).

    Parameters for the following functions:

    N

    The number of elements in the vectors.

    X

    Vector X.

    incX

    Stride within X. For example, if incX is 7, every 7th element is used.

    Y

    Vector Y.

    incY

    Stride within Y. For example, if incY is 7, every 7th element is used.

    returns

    See description above.

  100. def cblas_ssymm(Order: CBLAS_ORDER, Side: CBLAS_SIDE, Uplo: CBLAS_UPLO, M: CInt, N: CInt, alpha: CFloat, A: Ptr[CFloat], lda: CInt, B: Ptr[CFloat], ldb: CInt, beta: CFloat, C: Ptr[CFloat], ldc: CInt): Unit

    Permalink
  101. def cblas_ssymv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: CFloat, A: Ptr[CFloat], lda: CInt, X: Ptr[CFloat], incX: CInt, beta: CFloat, Y: Ptr[CFloat], incY: CInt): Unit

    Permalink
  102. def cblas_ssyr(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: CFloat, X: Ptr[CFloat], incX: CInt, A: Ptr[CFloat], lda: CInt): Unit

    Permalink
  103. def cblas_ssyr2(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: CFloat, X: Ptr[CFloat], incX: CInt, Y: Ptr[CFloat], incY: CInt, A: Ptr[CFloat], lda: CInt): Unit

    Permalink
  104. def cblas_ssyr2k(Order: CBLAS_ORDER, Uplo: CBLAS_UPLO, Trans: CBLAS_TRANSPOSE, N: CInt, K: CInt, alpha: CFloat, A: Ptr[CFloat], lda: CInt, B: Ptr[CFloat], ldb: CInt, beta: CFloat, C: Ptr[CFloat], ldc: CInt): Unit

    Permalink
  105. def cblas_ssyrk(Order: CBLAS_ORDER, Uplo: CBLAS_UPLO, Trans: CBLAS_TRANSPOSE, N: CInt, K: CInt, alpha: CFloat, A: Ptr[CFloat], lda: CInt, beta: CFloat, C: Ptr[CFloat], ldc: CInt): Unit

    Permalink
  106. def cblas_stbmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, K: CInt, A: Ptr[CFloat], lda: CInt, X: Ptr[CFloat], incX: CInt): Unit

    Permalink
  107. def cblas_stbsv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, K: CInt, A: Ptr[CFloat], lda: CInt, X: Ptr[CFloat], incX: CInt): Unit

    Permalink
  108. def cblas_stpmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, Ap: Ptr[CFloat], X: Ptr[CFloat], incX: CInt): Unit

    Permalink
  109. def cblas_stpsv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, Ap: Ptr[CFloat], X: Ptr[CFloat], incX: CInt): Unit

    Permalink
  110. def cblas_strmm(Order: CBLAS_ORDER, Side: CBLAS_SIDE, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, M: CInt, N: CInt, alpha: CFloat, A: Ptr[CFloat], lda: CInt, B: Ptr[CFloat], ldb: CInt): Unit

    Permalink
  111. def cblas_strmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, A: Ptr[CFloat], lda: CInt, X: Ptr[CFloat], incX: CInt): Unit

    Permalink
  112. def cblas_strsm(Order: CBLAS_ORDER, Side: CBLAS_SIDE, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, M: CInt, N: CInt, alpha: CFloat, A: Ptr[CFloat], lda: CInt, B: Ptr[CFloat], ldb: CInt): Unit

    Permalink
  113. def cblas_strsv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, A: Ptr[CFloat], lda: CInt, X: Ptr[CFloat], incX: CInt): Unit

    Permalink
  114. def cblas_xerbla(p: CInt, rout: CString, form: CString, varArgs: CVararg*): Unit

    Permalink

    The error handler for the LAPACK routines.

    The error handler for the LAPACK routines. It is called by an LAPACK routine if an input parameter has an invalid value. A message is printed and execution stops.

  115. def cblas_zaxpy(N: CInt, alpha: Ptr[CDoubleComplex], X: Ptr[CDoubleComplex], incX: CInt, Y: Ptr[CDoubleComplex], incY: CInt): Unit

    Permalink

    A constant times a vector plus a vector (double-precision complex).

  116. def cblas_zcopy(N: CInt, X: Ptr[CDoubleComplex], incX: CInt, Y: Ptr[CDoubleComplex], incY: CInt): Unit

    Permalink

    Copies a vector to another vector (double-precision complex).

  117. def cblas_zdotc_sub(N: CInt, X: Ptr[CDoubleComplex], incX: CInt, Y: Ptr[CDoubleComplex], incY: CInt, dotc: Ptr[CDoubleComplex]): Unit

    Permalink

    Dot product of the complex conjugate of a double-precision complex vector with a second double-precision complex vector.

    Dot product of the complex conjugate of a double-precision complex vector with a second double-precision complex vector.

    dotc

    - The result vector. Computes conjg(X) * Y.

  118. def cblas_zdotu_sub(N: CInt, X: Ptr[CDoubleComplex], incX: CInt, Y: Ptr[CDoubleComplex], incY: CInt, dotu: Ptr[CDoubleComplex]): Unit

    Permalink

    Dot product of two double-precision complex vectors.

    Dot product of two double-precision complex vectors.

    dotu

    - The result vector.

  119. def cblas_zdscal(N: CInt, alpha: CDouble, X: Ptr[CDoubleComplex], incX: CInt): Unit

    Permalink

    Multiplies each element of a vector by a constant (double-precision complex).

  120. def cblas_zgbmv(order: CBLAS_ORDER, TransA: CBLAS_TRANSPOSE, M: CInt, N: CInt, KL: CInt, KU: CInt, alpha: Ptr[CDoubleComplex], A: Ptr[CDoubleComplex], lda: CInt, X: Ptr[CDoubleComplex], incX: CInt, beta: Ptr[CDoubleComplex], Y: Ptr[CDoubleComplex], incY: CInt): Unit

    Permalink
  121. def cblas_zgemm(Order: CBLAS_ORDER, TransA: CBLAS_TRANSPOSE, TransB: CBLAS_TRANSPOSE, M: CInt, N: CInt, K: CInt, alpha: Ptr[CDoubleComplex], A: Ptr[CDoubleComplex], lda: CInt, B: Ptr[CDoubleComplex], ldb: CInt, beta: Ptr[CDoubleComplex], C: Ptr[CDoubleComplex], ldc: CInt): Unit

    Permalink
  122. def cblas_zgemv(order: CBLAS_ORDER, TransA: CBLAS_TRANSPOSE, M: CInt, N: CInt, alpha: Ptr[CDoubleComplex], A: Ptr[CDoubleComplex], lda: CInt, X: Ptr[CDoubleComplex], incX: CInt, beta: Ptr[CDoubleComplex], Y: Ptr[CDoubleComplex], incY: CInt): Unit

    Permalink
  123. def cblas_zgerc(order: CBLAS_ORDER, M: CInt, N: CInt, alpha: Ptr[CDoubleComplex], X: Ptr[CDoubleComplex], incX: CInt, Y: Ptr[CDoubleComplex], incY: CInt, A: Ptr[CDoubleComplex], lda: CInt): Unit

    Permalink
  124. def cblas_zgeru(order: CBLAS_ORDER, M: CInt, N: CInt, alpha: Ptr[CDoubleComplex], X: Ptr[CDoubleComplex], incX: CInt, Y: Ptr[CDoubleComplex], incY: CInt, A: Ptr[CDoubleComplex], lda: CInt): Unit

    Permalink
  125. def cblas_zhbmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, K: CInt, alpha: Ptr[CDoubleComplex], A: Ptr[CDoubleComplex], lda: CInt, X: Ptr[CDoubleComplex], incX: CInt, beta: Ptr[CDoubleComplex], Y: Ptr[CDoubleComplex], incY: CInt): Unit

    Permalink
  126. def cblas_zhemm(Order: CBLAS_ORDER, Side: CBLAS_SIDE, Uplo: CBLAS_UPLO, M: CInt, N: CInt, alpha: Ptr[CDoubleComplex], A: Ptr[CDoubleComplex], lda: CInt, B: Ptr[CDoubleComplex], ldb: CInt, beta: Ptr[CDoubleComplex], C: Ptr[CDoubleComplex], ldc: CInt): Unit

    Permalink
  127. def cblas_zhemv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: Ptr[CDoubleComplex], A: Ptr[CDoubleComplex], lda: CInt, X: Ptr[CDoubleComplex], incX: CInt, beta: Ptr[CDoubleComplex], Y: Ptr[CDoubleComplex], incY: CInt): Unit

    Permalink
  128. def cblas_zher(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: CDouble, X: Ptr[CDoubleComplex], incX: CInt, A: Ptr[CDoubleComplex], lda: CInt): Unit

    Permalink
  129. def cblas_zher2(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: Ptr[CDoubleComplex], X: Ptr[CDoubleComplex], incX: CInt, Y: Ptr[CDoubleComplex], incY: CInt, A: Ptr[CDoubleComplex], lda: CInt): Unit

    Permalink
  130. def cblas_zher2k(Order: CBLAS_ORDER, Uplo: CBLAS_UPLO, Trans: CBLAS_TRANSPOSE, N: CInt, K: CInt, alpha: Ptr[CDoubleComplex], A: Ptr[CDoubleComplex], lda: CInt, B: Ptr[CDoubleComplex], ldb: CInt, beta: CDouble, C: Ptr[CDoubleComplex], ldc: CInt): Unit

    Permalink
  131. def cblas_zherk(Order: CBLAS_ORDER, Uplo: CBLAS_UPLO, Trans: CBLAS_TRANSPOSE, N: CInt, K: CInt, alpha: CDouble, A: Ptr[CDoubleComplex], lda: CInt, beta: CDouble, C: Ptr[CDoubleComplex], ldc: CInt): Unit

    Permalink
  132. def cblas_zhpmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: Ptr[CDoubleComplex], Ap: Ptr[CDoubleComplex], X: Ptr[CDoubleComplex], incX: CInt, beta: Ptr[CDoubleComplex], Y: Ptr[CDoubleComplex], incY: CInt): Unit

    Permalink
  133. def cblas_zhpr(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: CDouble, X: Ptr[CDoubleComplex], incX: CInt, A: Ptr[CDoubleComplex]): Unit

    Permalink
  134. def cblas_zhpr2(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, N: CInt, alpha: Ptr[CDoubleComplex], X: Ptr[CDoubleComplex], incX: CInt, Y: Ptr[CDoubleComplex], incY: CInt, Ap: Ptr[CDoubleComplex]): Unit

    Permalink
  135. def cblas_zscal(N: CInt, alpha: Ptr[CDoubleComplex], X: Ptr[CDoubleComplex], incX: CInt): Unit

    Permalink

    Multiplies each element of a vector by a constant (double-precision complex).

  136. def cblas_zswap(N: CInt, X: Ptr[CDoubleComplex], incX: CInt, Y: Ptr[CDoubleComplex], incY: CInt): Unit

    Permalink

    Exchanges the elements of two vectors (double-precision complex).

  137. def cblas_zsymm(Order: CBLAS_ORDER, Side: CBLAS_SIDE, Uplo: CBLAS_UPLO, M: CInt, N: CInt, alpha: Ptr[CDoubleComplex], A: Ptr[CDoubleComplex], lda: CInt, B: Ptr[CDoubleComplex], ldb: CInt, beta: Ptr[CDoubleComplex], C: Ptr[CDoubleComplex], ldc: CInt): Unit

    Permalink
  138. def cblas_zsyr2k(Order: CBLAS_ORDER, Uplo: CBLAS_UPLO, Trans: CBLAS_TRANSPOSE, N: CInt, K: CInt, alpha: Ptr[CDoubleComplex], A: Ptr[CDoubleComplex], lda: CInt, B: Ptr[CDoubleComplex], ldb: CInt, beta: Ptr[CDoubleComplex], C: Ptr[CDoubleComplex], ldc: CInt): Unit

    Permalink
  139. def cblas_zsyrk(Order: CBLAS_ORDER, Uplo: CBLAS_UPLO, Trans: CBLAS_TRANSPOSE, N: CInt, K: CInt, alpha: Ptr[CDoubleComplex], A: Ptr[CDoubleComplex], lda: CInt, beta: Ptr[CDoubleComplex], C: Ptr[CDoubleComplex], ldc: CInt): Unit

    Permalink
  140. def cblas_ztbmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, K: CInt, A: Ptr[CDoubleComplex], lda: CInt, X: Ptr[CDoubleComplex], incX: CInt): Unit

    Permalink
  141. def cblas_ztbsv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, K: CInt, A: Ptr[CDoubleComplex], lda: CInt, X: Ptr[CDoubleComplex], incX: CInt): Unit

    Permalink
  142. def cblas_ztpmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, Ap: Ptr[CDoubleComplex], X: Ptr[CDoubleComplex], incX: CInt): Unit

    Permalink
  143. def cblas_ztpsv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, Ap: Ptr[CDoubleComplex], X: Ptr[CDoubleComplex], incX: CInt): Unit

    Permalink
  144. def cblas_ztrmm(Order: CBLAS_ORDER, Side: CBLAS_SIDE, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, M: CInt, N: CInt, alpha: Ptr[CDoubleComplex], A: Ptr[CDoubleComplex], lda: CInt, B: Ptr[CDoubleComplex], ldb: CInt): Unit

    Permalink
  145. def cblas_ztrmv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, A: Ptr[CDoubleComplex], lda: CInt, X: Ptr[CDoubleComplex], incX: CInt): Unit

    Permalink
  146. def cblas_ztrsm(Order: CBLAS_ORDER, Side: CBLAS_SIDE, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, M: CInt, N: CInt, alpha: Ptr[CDoubleComplex], A: Ptr[CDoubleComplex], lda: CInt, B: Ptr[CDoubleComplex], ldb: CInt): Unit

    Permalink
  147. def cblas_ztrsv(order: CBLAS_ORDER, Uplo: CBLAS_UPLO, TransA: CBLAS_TRANSPOSE, Diag: CBLAS_DIAG, N: CInt, A: Ptr[CDoubleComplex], lda: CInt, X: Ptr[CDoubleComplex], incX: CInt): Unit

    Permalink
  148. def clone(): AnyRef

    Permalink
    Attributes
    protected[java.lang]
    Definition Classes
    AnyRef
    Annotations
    @throws( ... )
  149. final def eq(arg0: AnyRef): Boolean

    Permalink
    Definition Classes
    AnyRef
  150. def equals(arg0: Any): Boolean

    Permalink
    Definition Classes
    AnyRef → Any
  151. def finalize(): Unit

    Permalink
    Attributes
    protected[java.lang]
    Definition Classes
    AnyRef
    Annotations
    @throws( classOf[java.lang.Throwable] )
  152. final def getClass(): Class[_]

    Permalink
    Definition Classes
    AnyRef → Any
  153. def hashCode(): Int

    Permalink
    Definition Classes
    AnyRef → Any
  154. final def isInstanceOf[T0]: Boolean

    Permalink
    Definition Classes
    Any
  155. final def ne(arg0: AnyRef): Boolean

    Permalink
    Definition Classes
    AnyRef
  156. final def notify(): Unit

    Permalink
    Definition Classes
    AnyRef
  157. final def notifyAll(): Unit

    Permalink
    Definition Classes
    AnyRef
  158. final def synchronized[T0](arg0: ⇒ T0): T0

    Permalink
    Definition Classes
    AnyRef
  159. def toString(): String

    Permalink
    Definition Classes
    AnyRef → Any
  160. final def wait(): Unit

    Permalink
    Definition Classes
    AnyRef
    Annotations
    @throws( ... )
  161. final def wait(arg0: Long, arg1: Int): Unit

    Permalink
    Definition Classes
    AnyRef
    Annotations
    @throws( ... )
  162. final def wait(arg0: Long): Unit

    Permalink
    Definition Classes
    AnyRef
    Annotations
    @throws( ... )

Inherited from AnyRef

Inherited from Any

Ungrouped