.hx.
Computes conjugate transpose matrix-matrix product.
Operator Return Value
Matrix containing the product of AH and B. (Output)
Required Operands
A — Left operand matrix. This is an array of rank 2 or 3. It may be real, double, complex, double complex, or one of the computational sparse matrix derived types, c_hbc_sparse or z_hbc_sparse. (Input)
Note that A and B cannot both be ?_hbc_sparse.
B — Right operand matrix or vector. This is an array of rank 1, 2, or 3. It may be real, double, complex, double complex, or one of the computational sparse matrix derived types, c_hbc_sparse or z_hbc_sparse. (Input)
Note that A and B cannot both be ?_hbc_sparse.
If A has rank three, B must have rank three.
If B has rank three, A must have rank three.
FORTRAN 90 Interface
A .hx. B
Description
Computes the product of the conjugate transpose of matrix A and matrix or vector B. The results are in a precision and data type that ascends to the most accurate or complex operand.
Rank three operation is defined as follows:
do i = 1, min(size(A,3), size(B,3))
X(:,:,i) = A(:,:,i) .hx. B(:,:,i)
end do
.hx. can be used with either dense or sparse matrices. It is MPI capable for dense matrices only.
Examples
Dense Matrix Example (operator_ex32.f90)
use linear_operators
implicit none
! This is the equivalent of Example 4 (using operators) for LIN_EIG_GEN.
integer, parameter :: n=17
real(kind(1d0)), parameter :: one=1d0
real(kind(1d0)), dimension(n,n) :: A, C
real(kind(1d0)) variation(n), eta
complex(kind(1d0)), dimension(n,n) :: U, V, e(n), d(n)
! Generate a random matrix.
A = rand(A)
! Compute the eigenvalues, left- and right- eigenvectors.
D = EIG(A, W=V); E = EIG(.t.A, W=U)
! Compute condition numbers and variations of eigenvalues.
variation = norm(A)/abs(diagonals( U .hx. V))
! Now perturb the data in the matrix by the relative factors
! eta=sqrt(epsilon) and solve for values again. Check the
! differences compared to the estimates. They should not exceed
! the bounds.
eta = sqrt(epsilon(one))
C = A + eta*(2*rand(A)-1)*A
D = EIG(C)
! Looking at the differences of absolute values accounts for
! switching signs on the imaginary parts.
if (count(abs(d)-abs(e) > eta*variation) == 0) then
write (*,*) 'Example 4 for LIN_EIG_GEN (operators) is correct.'
end if
end
use wrcrn_int
use linear_operators
type (c_sparse) S
type (c_hbc_sparse) H
integer, parameter :: N=3
complex (kind(1.e0)) x(N,N), y(N,N), A(N,N)
real (kind(1.e0)) err
S = c_entry (1, 1, (2.0, 1.0) )
S = c_entry (1, 3, (1.0, 3.0))
S = c_entry (2, 2, (4.0, -1.0))
S = c_entry (3, 3, (6.0, 2.0))
H = S ! sparse
X = H ! dense equivalent of H
A= rand(A)
Y = H .hx. A
call wrcrn ( 'H', X)
call wrcrn ( 'A', a)
call wrcrn ( 'H .hx. A ', y)
! Check the results.
err = norm(y - (X .hx. A))
if (err <= sqrt(epsilon(one))) then
write (*,*) 'Sparse example for .hx. operator is correct.'
end if
end
Output
H
1 2 3
1 ( 2.000, 1.000) ( 0.000, 0.000) ( 1.000, 3.000)
2 ( 0.000, 0.000) ( 4.000,-1.000) ( 0.000, 0.000)
3 ( 0.000, 0.000) ( 0.000, 0.000) ( 6.000, 2.000)
A
1 2 3
1 ( 0.6278, 0.8475) ( 0.8007, 0.4179) ( 0.4512, 0.2601)
2 ( 0.1249, 0.4675) ( 0.7957, 0.1609) ( 0.4228, 0.0507)
3 ( 0.4608, 0.0891) ( 0.3181, 0.9180) ( 0.9961, 0.1939)
H .hx. A
1 2 3
1 ( 2.103, 1.067) ( 2.019, 0.035) ( 1.163, 0.069)
2 ( 0.032, 1.995) ( 3.022, 1.439) ( 1.640, 0.626)
3 ( 6.113,-1.423) ( 5.799, 2.888) ( 7.596,-1.922)
Sparse example for .hx. operator is correct.
use linear_operators
use mpi_setup_int
integer, parameter :: N=32, nr=4
complex (kind(1.e0)) A(N,N,nr), B(N,N,nr), Y(N,N,nr)
! Setup for MPI
mp_nprocs = mp_setup()
if (mp_rank == 0) then
A = rand(A)
B = rand(B)
end if
Y = A .hx. B
mp_nprocs = mp_setup ('Final')
end