.tx.

   more...

   more...
Computes transpose matrix-matrix or transpose matrix-vector product.
Operator Return Value
Matrix containing the product of AT 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, ?_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, ?_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 .tx. B
Description
Computes the product of the 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) .tx. B(:,:,i)
end do
.tx. can be used with either dense or sparse matrices. It is MPI capable for dense matrices only.
Examples
Dense Matrix Example (operator_ex05.f90)
 
use linear_operators
implicit none
! This is the equivalent of Example 1 for LIN_SOL_SELF using operators
! and functions.
integer, parameter :: m=64, n=32
real(kind(1e0)) :: one=1.0e0, err
real(kind(1e0)) A(n,n), b(n,n), C(m,n), d(m,n), x(n,n)
! Generate two rectangular random matrices.
C = rand(C); d=rand(d)
! Form the normal equations for the rectangular system.
A = C .tx. C; b = C .tx. d
! Compute the solution for Ax = b, A is symmetric.
x = A .ix. b
! Check the results.
err = norm(b - (A .x. x))/(norm(A)+norm(b))
if (err <= sqrt(epsilon(one))) then
write (*,*) 'Example 1 for LIN_SOL_SELF (operators) is correct.'
end if
end
Sparse Matrix Example
 
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), B(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
B = rand(B)
Y = H .tx. B
call wrrrn ( 'H', X)
call wrrrn ( 'B', b)
call wrrrn ( 'H .tx. B ', y)
 
! Check the results.
err = norm(y - (X .tx. B))
if (err <= sqrt(epsilon(one))) then
write (*,*) 'Sparse example for .tx. operator is correct.'
end if
 
end
Output
 
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
B
1 2 3
1 0.8711 0.4467 0.4743
2 0.8315 0.7257 0.4518
3 0.6839 0.0561 0.6972
H .tx. B
1 2 3
1 1.742 0.893 0.949
2 3.326 2.903 1.807
3 4.975 0.784 4.657
Sparse example for .tx. operator is correct.
Parallel Example (parallel_ex05.f90)
 
use linear_operators
use mpi_setup_int
 
implicit none
 
! This is the equivalent of Parallel Example 5 using box data types,
! operators and functions.
 
integer, parameter :: m=64, n=32, nr=4
real(kind(1e0)) :: one=1e0, err(nr)
real(kind(1e0)), dimension(n,n,nr) :: A, b, x
real(kind(1e0)), dimension(m,n,nr) :: C, d
 
! Setup for MPI.
mp_nprocs = mp_setup()
 
! Generate two rectangular random matrices, only
! at the root node.
if (mp_rank == 0) then
C = rand(C); d=rand(d)
endif
 
! Form the normal equations for the rectangular system.
A = C .tx. C; b = C .tx. d
 
! Compute the solution for Ax = b.
x = A .ix. b
 
! Check the results.
err = norm(b - (A .x. x))/(norm(A)+norm(b))
if (ALL(err <= sqrt(epsilon(one))) .AND. MP_RANK == 0) &
write (*,*) 'Parallel Example 5 is correct.'
! See to any error messages and quit MPI.
mp_nprocs = mp_setup('Final')
 
end