.h.

Computes the conjugate transpose of a matrix.

Operator Return Value

Matrix containing the conjugate transpose of A.   (Output)

Required Operand

A — Matrix for which the conjugate transpose is to be computed. 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)

FORTRAN 90 Interface

.h. A

Description

Computes the conjugate transpose of  matrix A. The operation may be read adjoint, and the results are the mathematical objects in a precision and data type that matches the operand.  Since this is a unary operation, it has higher Fortran 90 precedence than any other intrinsic unary array operation. 

.h. can be used with either dense or sparse matrices.

Examples

Dense Matrix Example  (operator_ex34.f90)

 

      use linear_operators

 

      implicit none

 

! This is the equivalent of Example 2 (using operators) for LIN_GEIG_GEN.

 

      integer, parameter :: n=32

      real(kind(1d0)), parameter :: one=1d0, zero=0d0

      real(kind(1d0)) err, alpha(n)

      complex(kind(1d0)), dimension(n,n) :: A, B, C, D, V

 

 

! Generate random matrices for both A and B.

      C = rand(C); D = rand(D)

      A = C + .h.C; B = D .hx. D; B = (B + .h.B)/2

 

      ALPHA = EIG(A, B=B, W=V)

 

! Check that residuals are small.  Use a real array for  alpha 

! since the eigenvalues are known to be real.

      err= norm((A .x. V) - (B .x. V .x. diag(alpha)),1)/&

           (norm(A,1)+norm(B,1)*norm(alpha,1))

      if (err <= sqrt(epsilon(one))) then

         write (*,*) 'Example 2 for LIN_GEIG_GEN (operators) is correct.'

      end if

 

      end 

Sparse Matrix Example

 

 use wrcrn_int

 use linear_operators

 

 type (c_sparse) S

 type (c_hbc_sparse) H, HT

 integer, parameter :: N=3

 complex (kind(1.e0)) X(3,3), XT(3,3)

 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

 HT = .h. H

 XT = HT ! dense equivalent of HT

  call wrcrn ( 'H', X)

  call wrcrn ( 'H Conjugate Transpose', XT)

 

! Check the results.

     err =  norm(XT - (.h. X))

      if (err <= sqrt(epsilon(one))) then

         write (*,*) 'Sparse example for .h. 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)

                 H Conjugate Transpose

                  1                2                3

 1  ( 2.000,-1.000)  ( 0.000, 0.000)  ( 0.000, 0.000)

 2  ( 0.000, 0.000)  ( 4.000, 1.000)  ( 0.000, 0.000)

 3  ( 1.000,-3.000)  ( 0.000, 0.000)  ( 6.000,-2.000)

 Sparse example for .h. operator is correct.


Visual Numerics, Inc.
Visual Numerics - Developers of IMSL and PV-WAVE
http://www.vni.com/
PHONE: 713.784.3131
FAX:713.781.9260