.xh.

   more...

   more...
Computes a matrix-conjugate transpose matrix product.
Operator Return Value
Matrix containing the product of A and BH. (Output)
Required Operands
A — Left 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.
B — Right 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.
If A has rank three, B must have rank three.
If B has rank three,
A must have rank three.
FORTRAN 90 Interface
A .xh. B
Description
Computes the product of matrix or vector A and the conjugate transpose of matrix 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) .xh. B(:,:,i)
end do
.xh. can be used with either dense or sparse matrices. It is MPI capable for dense matrices only.
Examples
Dense Matrix Example
 
use wrcrn_int
use linear_operators
integer, parameter :: N=3
complex (kind(1.e0)) A(N,N), B(N,N), Y(N,N)
 
A = rand(A)
B = rand(B)
Y = A .xh. B
call wrcrn ( 'A', a)
call wrcrn ( 'B', b)
call wrcrn ( 'A .xh. B ', y)
end
Output
A
1 2 3
1 ( 0.8071, 0.0054) ( 0.5617, 0.2508) ( 0.0223, 0.5555)
2 ( 0.9380, 0.5181) ( 0.8895, 0.9512) ( 0.7951, 0.6010)
3 ( 0.8349, 0.7291) ( 0.4162, 0.5255) ( 0.7388, 0.0309)
B
1 2 3
1 ( 0.5342, 0.2246) ( 0.9045, 0.0550) ( 0.4576, 0.3173)
2 ( 0.5531, 0.3362) ( 0.0757, 0.3970) ( 0.6807, 0.8625)
3 ( 0.3553, 0.9157) ( 0.0951, 0.7807) ( 0.4853, 0.0617)
A .xh. B
1 2 3
1 ( 1.141, 0.265) ( 1.085,-0.113) ( 0.586,-0.884)
2 ( 2.029, 0.900) ( 2.198,-0.587) ( 2.058,-1.036)
3 ( 1.363, 0.434) ( 1.477,-0.619) ( 1.775,-0.811)
Spare Matrix Example
 
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 = A .xh. H
call wrcrn ( 'A', a)
call wrcrn ( 'H', X)
call wrcrn ( 'A .xh. H ', y)
 
! Check the results.
err = norm(y - (A .xh. X))
if (err <= sqrt(epsilon(one))) then
write (*,*) 'Sparse example for .xh. operator is correct.'
end if
 
end
Output
A
1 2 3
1 ( 0.8526, 0.3532) ( 0.1822, 0.3938) ( 0.8008, 0.1308)
2 ( 0.5599, 0.8914) ( 0.7541, 0.5163) ( 0.8713, 0.9580)
3 ( 0.9947, 0.2735) ( 0.6237, 0.2137) ( 0.3802, 0.8903)
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 .xh. H
1 2 3
1 ( 3.252,-2.418) ( 0.335, 1.757) ( 5.066,-0.817)
2 ( 5.757,-0.433) ( 2.500, 2.819) ( 7.144, 4.005)
3 ( 5.314,-0.698) ( 2.281, 1.478) ( 4.062, 4.581)
Sparse example for .xh. operator is correct.
Parallel Example
 
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 .xh. B
mp_nprocs = mp_setup ('Final')
 
end