.xt.
Computes matrix-transpose matrix product.
Operator Return Value
Matrix containing the product of A and BT. (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, ?_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, ?_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 .xt. B
Description
Computes the product of matrix or vector A and the 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) .xt. B(:,:,i)
end do
.xt. can be used with either dense or sparse matrices. It is MPI capable for dense matrices only.
Examples
Dense Matrix Example (operator_ex14.f90)
use linear_operators
implicit none
!
integer, parameter :: n=32
real(kind(1d0)) :: one=1d0, zero=0d0
real(kind(1d0)) A(n,n), P(n,n), Q(n,n), &
S_D(n), U_D(n,n), V_D(n,n)
! Generate a random matrix.
A = rand(A)
! Compute the singular value decomposition.
S_D = SVD(A, U=U_D, V=V_D)
! Compute the (left) orthogonal factor.
P = U_D .xt. V_D
! Compute the (right) self-adjoint factor.
Q = V_D .x. diag(S_D) .xt. V_D
! Check the results.
if (norm( EYE(n) - (P .xt. P)) &
<= sqrt(epsilon(one))) then
if (norm(A - (P .x. Q))/norm(A) &
<= sqrt(epsilon(one))) then
write (*,*) 'Example 2 for LIN_SOL_SVD (operators) is correct.'
end if
end if
end
use wrrrn_int
use linear_operators
type (s_sparse) S
type (s_hbc_sparse) H
integer, parameter :: N=3
real (kind(1.e0)) x(N,N), y(N,N), a(N,N)
real (kind(1.e0)) err
S = s_entry (1, 1, 2.0)
S = s_entry (1, 3, 1.0)
S = s_entry (2, 2, 4.0)
S = s_entry (3, 3, 6.0)
H = S ! sparse
X = H ! dense equivalent of H
A = rand(A)
Y = A .xt. H
call wrrrn ( 'A', A)
call wrrrn ( 'H', X)
call wrrrn ( 'A .xt. H', y)
! Check the results.
err = norm(y - (A .xt. X))
if (err <= sqrt(epsilon(one))) then
write (*,*) 'Sparse example for .xt. operator is correct.'
end if
end
Output
A
1 2 3
1 0.5423 0.2380 0.9250
2 0.0844 0.1323 0.1937
3 0.4146 0.3135 0.7757
H
1 2 3
1 2.000 0.000 1.000
2 0.000 4.000 0.000
3 0.000 0.000 6.000
A .xt. H
1 2 3
1 2.010 0.952 5.550
2 0.363 0.529 1.162
3 1.605 1.254 4.654
Sparse example for .xt. operator is correct.
Parallel Example (parallel_ex15.f90)
A “Polar Decomposition” of several matrices are computed. The box data type and the SVD() function are used. Orthogonality and small residuals are checked to verify that the results are correct.
use linear_operators
use mpi_setup_int
implicit none
! This is the equivalent of Parallel Example 15 using operators and,
! functions for a polar decomposition.
integer, parameter :: n=33, nr=3
real(kind(1d0)) :: one=1d0, zero=0d0
real(kind(1d0)),dimension(n,n,nr) :: A, P, Q, &
S_D(n,nr), U_D, V_D
real(kind(1d0)) TEMP1(nr), TEMP2(nr)
! Setup for MPI:
mp_nprocs = mp_setup()
! Generate a random matrix.
if(mp_rank == 0) A = rand(A)
! Compute the singular value decomposition.
S_D = SVD(A, U=U_D, V=V_D)
! Compute the (left) orthogonal factor.
P = U_D .xt. V_D
! Compute the (right) self-adjoint factor.
Q = V_D .x. diag(S_D) .xt. V_D
! Check the results for orthogonality and
! small residuals.
TEMP1 = NORM(spread(EYE(n),3,nr) - (p .xt. p))
TEMP2 = NORM(A -(P .X. Q)) / NORM(A)
if (ALL(TEMP1 <= sqrt(epsilon(one))) .and. &
ALL(TEMP2 <= sqrt(epsilon(one)))) then
if(mp_rank == 0)&
write (*,*) 'Parallel Example 15 is correct.'
end if
! See to any error messages and exit MPI.
mp_nprocs = mp_setup('Final')
end