program main c*********************************************************************72 c cc blas3_d_test() tests blas3_d(). c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 06 April 2014 c c Author: c c John Burkardt c implicit none call timestamp ( ) write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'blas3_d_test():' write ( *, '(a)' ) ' FORTRAN77 version' write ( *, '(a)' ) ' Test blas3_d().' call dgemm_test ( ) call dtrmm_test ( ) call dtrsm_test ( ) c c Terminate. c write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'blas3_d_test():' write ( *, '(a)' ) ' Normal end of execution.' write ( *, '(a)' ) ' ' call timestamp ( ) stop end subroutine dgemm_test ( ) c********************************************************************72 c cc DGEMM_TEST tests DGEMM. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 10 February 2014 c c Author: c c John Burkardt c implicit none double precision a(10*10) double precision alpha double precision b(10*10) double precision beta double precision c(10*10) integer k integer lda integer ldb integer ldc integer m integer n character transa character transb character transc write ( *, '(a)' ) '' write ( *, '(a)' ) 'DGEMM_TEST' write ( *, '(a)' ) ' DGEMM carries out matrix multiplications' write ( *, '(a)' ) ' for double precision real matrices.' write ( *, '(a)' ) '' write ( *, '(a)' ) ' 1: C = alpha * A * B + beta * C;' write ( *, '(a)' ) ' 2: C = alpha * A'' * B + beta * C;' write ( *, '(a)' ) ' 3: C = alpha * A * B'' + beta * C;' write ( *, '(a)' ) ' 4: C = alpha * A'' * B'' + beta * C;' write ( *, '(a)' ) '' write ( *, '(a)' ) & ' We carry out all four calculations, but in each case,' write ( *, '(a)' ) & ' we choose our input matrices so that we get the same result.' c c C = alpha * A * B + beta * C. c transa = 'N' transb = 'N' transc = 'N' m = 4 n = 5 k = 3 alpha = 2.0D+00 lda = m call r8mat_test ( transa, lda, m, k, a ) ldb = k call r8mat_test ( transb, ldb, k, n, b ) beta = 3.0D+00 ldc = m call r8mat_test ( transc, ldc, m, n, c ) call dgemm ( transa, transb, m, n, k, alpha, a, lda, b, ldb, & beta, c, ldc ) call r8mat_print ( m, n, c, & ' C = alpha * A * B + beta * C:' ) c c C = alpha * A' * B + beta * C. c transa = 'T' transb = 'N' transc = 'N' m = 4 n = 5 k = 3 alpha = 2.0D+00 lda = k call r8mat_test ( transa, lda, m, k, a ) ldb = k call r8mat_test ( transb, ldb, k, n, b ) beta = 3.0D+00 ldc = m call r8mat_test ( transc, ldc, m, n, c ) call dgemm ( transa, transb, m, n, k, alpha, a, lda, b, ldb, & beta, c, ldc ) call r8mat_print ( m, n, c, & ' C = alpha * A'' * B + beta * C:' ) c c C = alpha * A * B' + beta * C. c transa = 'N' transb = 'T' transc = 'N' m = 4 n = 5 k = 3 alpha = 2.0D+00 lda = m call r8mat_test ( transa, lda, m, k, a ) ldb = n call r8mat_test ( transb, ldb, k, n, b ) beta = 3.0D+00 ldc = m call r8mat_test ( transc, ldc, m, n, c ) call dgemm ( transa, transb, m, n, k, alpha, a, lda, b, ldb, & beta, c, ldc ) call r8mat_print ( m, n, c, & ' C = alpha * A * B'' + beta * C:' ) c c C = alpha * A' * B' + beta * C. c transa = 'T' transb = 'T' transc = 'N' m = 4 n = 5 k = 3 alpha = 2.0D+00 lda = k call r8mat_test ( transa, lda, m, k, a ) ldb = n call r8mat_test ( transb, ldb, k, n, b ) beta = 3.0D+00 ldc = m call r8mat_test ( transc, ldc, m, n, c ) call dgemm ( transa, transb, m, n, k, alpha, a, lda, b, ldb, & beta, c, ldc ) call r8mat_print ( m, n, c, & ' C = alpha * A'' * B'' + beta * C:' ) return end subroutine dtrmm_test ( ) c*********************************************************************72 c cc DTRMM_TEST tests DTRMM. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 05 April 2014 c c Author: c c John Burkardt c implicit none integer m parameter ( m = 4 ) integer n parameter ( n = 5 ) integer lda parameter ( lda = m ) integer ldb parameter ( ldb = m ) double precision a(lda,m) double precision alpha double precision b(ldb,n) character diag integer i integer j character side character transa character transb character uplo write ( *, '(a)' ) '' write ( *, '(a)' ) 'DTRMM_TEST' write ( *, '(a)' ) & ' DTRMM multiplies a triangular matrix A and a' write ( *, '(a)' ) ' rectangular matrix B' write ( *, '(a)' ) '' write ( *, '(a)' ) ' 1: B = alpha * A * B;' write ( *, '(a)' ) ' 2: B = alpha * A'' * B;' c c B = alpha * A * B. c side = 'L' uplo = 'U' transa = 'N' diag = 'N' alpha = 2.0D+00 do j = 1, m do i = 1, j a(i,j) = dble ( i + j ) end do do i = j + 1, m a(i,j) = 0.0D+00 end do end do transb = 'N' call r8mat_test ( transb, ldb, m, n, b ) call dtrmm ( side, uplo, transa, diag, m, n, alpha, a, lda, b, & ldb ) call r8mat_print ( m, n, b, ' B = alpha * A * B:' ); c c B = alpha * A' * B. c side = 'L' uplo = 'U' transa = 'T' diag = 'N' alpha = 2.0D+00 do j = 1, m do i = 1, j a(i,j) = dble ( i + j ) end do do i = j + 1, m a(i,j) = 0.0D+00 end do end do transb = 'N' call r8mat_test ( transb, ldb, m, n, b ) call dtrmm ( side, uplo, transa, diag, m, n, alpha, a, lda, b, & ldb ) call r8mat_print ( m, n, b, ' B = alpha * A * B:' ); return end subroutine dtrsm_test ( ) c*********************************************************************72 c cc DTRSM_TEST tests DTRSM. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 06 April 2014 c c Author: c c John Burkardt c implicit none integer m parameter ( m = 4 ) integer n parameter ( n = 5 ) double precision amm(m,m) double precision ann(n,n) double precision alpha double precision b(m,n) character diag integer i integer j integer lda integer ldb character side character transa character transb character uplo write ( *, '(a)' ) '' write ( *, '(a)' ) 'DTRSM_TEST' write ( *, '(a)' ) & ' DTRSM solves a linear system involving a triangular' write ( *, '(a)' ) ' matrix A and a rectangular matrix B.' write ( *, '(a)' ) '' write ( *, '(a)' ) ' 1: Solve A * X = alpha * B;' write ( *, '(a)' ) ' 2: Solve A'' * X = alpha * B;' write ( *, '(a)' ) ' 3: Solve X * A = alpha * B;' write ( *, '(a)' ) ' 4: Solve X * A'' = alpha * B;' c c Solve A * X = alpha * B. c side = 'L' uplo = 'U' transa = 'N' diag = 'N' alpha = 2.0D+00 lda = m ldb = m do j = 1, m do i = 1, j amm(i,j) = dble ( i + j ) end do do i = j + 1, m amm(i,j) = 0.0D+00 end do end do transb = 'N' call r8mat_test ( transb, ldb, m, n, b ) call dtrsm ( side, uplo, transa, diag, m, n, alpha, amm, lda, b, & ldb ) call r8mat_print ( m, n, b, ' X = inv ( A ) * alpha * B:' ); c c Solve A' * X = alpha * B. c side = 'L' uplo = 'U' transa = 'T' diag = 'N' alpha = 2.0D+00 lda = m ldb = m do j = 1, m do i = 1, j amm(i,j) = dble ( i + j ) end do do i = j + 1, m amm(i,j) = 0.0D+00 end do end do transb = 'N' call r8mat_test ( transb, ldb, m, n, b ) call dtrsm ( side, uplo, transa, diag, m, n, alpha, amm, lda, b, & ldb ) call r8mat_print ( m, n, b, ' X = inv ( A'' ) * alpha * B:' ); c c Solve X * A = alpha * B. c side = 'R' uplo = 'U' transa = 'N' diag = 'N' alpha = 2.0D+00 lda = n ldb = m do j = 1, n do i = 1, j ann(i,j) = dble ( i + j ) end do do i = j + 1, n ann(i,j) = 0.0D+00 end do end do transb = 'N' call r8mat_test ( transb, ldb, m, n, b ) call dtrsm ( side, uplo, transa, diag, m, n, alpha, ann, lda, b, & ldb ) call r8mat_print ( m, n, b, ' X = alpha * B * inv ( A ):' ); c c Solve X * A'' = alpha * B. c side = 'R' uplo = 'U' transa = 'T' diag = 'N' alpha = 2.0D+00 lda = n ldb = m do j = 1, n do i = 1, j ann(i,j) = dble ( i + j ) end do do i = j + 1, n ann(i,j) = 0.0D+00 end do end do transb = 'N' call r8mat_test ( transb, ldb, m, n, b ) call dtrsm ( side, uplo, transa, diag, m, n, alpha, ann, lda, b, & ldb ) call r8mat_print ( m, n, b, ' X = alpha * B * inv ( A'' ):' ); return end function lsame ( ca, cb ) c*********************************************************************72 c cc lsame() returns TRUE if CA is the same letter as CB regardless of case. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 09 February 2014 c c Author: c c Original FORTRAN77 version by Jack Dongarra. c This version by John Burkardt. c c Reference: c c Jack Dongarra, Jim Bunch, Cleve Moler, Pete Stewart, c LINPACK User's Guide, c SIAM, 1979, c ISBN13: 978-0-898711-72-1, c LC: QA214.L56. c c Charles Lawson, Richard Hanson, David Kincaid, Fred Krogh, c Basic Linear Algebra Subprograms for FORTRAN usage, c ACM Transactions on Mathematical Software, c Volume 5, Number 3, pages 308-323, 1979. c c Parameters: c c Input, character CA, CB, the character to compare. c c Output, logical LSAME, is TRUE if the characters are equal, c disregarding case. c implicit none character ca character cb intrinsic ichar integer inta integer intb logical lsame integer zcode c c Test if the characters are equal. c lsame = ca .eq. cb if ( lsame ) then return end if c c Now test for equivalence if both characters are alphabetic. c zcode = ichar ( 'Z' ) c c Use 'Z' rather than 'A' so that ASCII can be detected on Prime c machines, on which ICHAR returns a value with bit 8 set. c ICHAR('A') on Prime machines returns 193 which is the same as c ICHAR('A') on an EBCDIC machine. c inta = ichar ( ca ) intb = ichar ( cb ) c c ASCII is assumed. c ZCODE is the ASCII code of either lower or upper case 'Z'. c if ( zcode .eq. 90 .or. zcode .eq. 122 ) then if ( inta .ge. 97 .and. inta .le. 122 ) then inta = inta - 32 end if if ( intb .ge. 97 .and. intb .le. 122 ) then intb = intb - 32 end if c c EBCDIC is assumed. c ZCODE is the EBCDIC code of either lower or upper case 'Z'. c else if ( zcode .eq. 233 .or. zcode .eq. 169 ) then if ( inta .ge. 129 .and. inta .le. 137 .or. & inta .ge. 145 .and. inta .le. 153 .or. & inta .ge. 162 .and. inta .le. 169 ) then inta = inta + 64 end if if ( intb .ge. 129 .and. intb .le. 137 .or. & intb .ge. 145 .and. intb .le. 153 .or. & intb .ge. 162 .and. intb .le. 169 ) then intb = intb + 64 end if c c ASCII is assumed, on Prime machines. c ZCODE is the ASCII code plus 128 of either lower or upper case 'Z'. c else if ( zcode .eq. 218 .or. zcode .eq. 250 ) then if ( inta .ge. 225 .and. inta .le. 250 ) then inta = inta - 32 end if if ( intb .ge. 225 .and. intb .le. 250 ) then intb = intb - 32 end if end if lsame = inta .eq. intb return end subroutine r8mat_print ( m, n, a, title ) c*********************************************************************72 c cc r8mat_print() prints an R8MAT. c c Discussion: c c An R8MAT is an array of R8's. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 20 May 2004 c c Author: c c John Burkardt c c Parameters: c c Input, integer M, the number of rows in A. c c Input, integer N, the number of columns in A. c c Input, double precision A(M,N), the matrix. c c Input, character ( len = * ) TITLE, a title. c implicit none integer m integer n double precision a(m,n) character ( len = * ) title call r8mat_print_some ( m, n, a, 1, 1, m, n, title ) return end subroutine r8mat_print_some ( m, n, a, ilo, jlo, ihi, jhi, & title ) c*********************************************************************72 c cc r8mat_print_some() prints some of an R8MAT. c c Discussion: c c An R8MAT is an array of R8's. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 25 January 2007 c c Author: c c John Burkardt c c Parameters: c c Input, integer M, N, the number of rows and columns. c c Input, double precision A(M,N), an M by N matrix to be printed. c c Input, integer ILO, JLO, the first row and column to print. c c Input, integer IHI, JHI, the last row and column to print. c c Input, character ( len = * ) TITLE, a title. c implicit none integer incx parameter ( incx = 5 ) integer m integer n double precision a(m,n) character * ( 14 ) ctemp(incx) integer i integer i2hi integer i2lo integer ihi integer ilo integer inc integer j integer j2 integer j2hi integer j2lo integer jhi integer jlo character * ( * ) title write ( *, '(a)' ) ' ' write ( *, '(a)' ) trim ( title ) if ( m .le. 0 .or. n .le. 0 ) then write ( *, '(a)' ) ' ' write ( *, '(a)' ) ' (None)' return end if do j2lo = max ( jlo, 1 ), min ( jhi, n ), incx j2hi = j2lo + incx - 1 j2hi = min ( j2hi, n ) j2hi = min ( j2hi, jhi ) inc = j2hi + 1 - j2lo write ( *, '(a)' ) ' ' do j = j2lo, j2hi j2 = j + 1 - j2lo write ( ctemp(j2), '(i7,7x)') j end do write ( *, '('' Col '',5a14)' ) ( ctemp(j), j = 1, inc ) write ( *, '(a)' ) ' Row' write ( *, '(a)' ) ' ' i2lo = max ( ilo, 1 ) i2hi = min ( ihi, m ) do i = i2lo, i2hi do j2 = 1, inc j = j2lo - 1 + j2 write ( ctemp(j2), '(g14.6)' ) a(i,j) end do write ( *, '(i5,a,5a14)' ) i, ':', ( ctemp(j), j = 1, inc ) end do end do return end subroutine r8mat_test ( trans, lda, m, n, a ) c*********************************************************************72 c cc r8mat_test() sets up a test matrix. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 10 February 2014 c c Author: c c John Burkardt. c c Parameters: c c Input, character * ( 1 ) TRANS, indicates whether matrix is to be c transposed. c 'N', no transpose. c 'T', transpose the matrix. c c Input, integer LDA, the leading dimension of the matrix. c c Input, integer M, N, the number of rows and columns of c the matrix. c c Output, double precision A(LDA,*), the matrix. c if TRANS is 'N', then the matrix is stored in LDA*N entries, c as an M x N matrix; c if TRANS is 'T', then the matrix is stored in LDA*M entries, c as an N x M matrix. c implicit none integer lda double precision a(lda,*) integer i integer j integer m integer n character * ( 1 ) trans if ( trans .eq. 'N' ) then do j = 1, n do i = 1, m a(i,j) = dble ( 10 * i + j ) end do end do else do j = 1, n do i = 1, m a(j,i) = dble ( 10 * i + j ) end do end do end if return end subroutine timestamp ( ) c*********************************************************************72 c cc timestamp() prints out the current YMDHMS date as a timestamp. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 16 September 2005 c c Author: c c John Burkardt c c Parameters: c c None c implicit none character ( len = 8 ) date character ( len = 10 ) time call date_and_time ( date, time ) write ( *, '(a8,2x,a10)' ) date, time return end subroutine xerbla ( srname, info ) c*********************************************************************72 c cc xerbla() is an error handler. c c Discussion: c c XERBLA is called if an input parameter has an invalid value. c A message is printed and execution stops. c c Licensing: c c This code is distributed under the MIT license. c c Modified: c c 09 February 2014 c c Author: c c Original FORTRAN77 version by Jack Dongarra. c This version by John Burkardt. c c Reference: c c Jack Dongarra, Jim Bunch, Cleve Moler, Pete Stewart, c LINPACK User's Guide, c SIAM, 1979, c ISBN13: 978-0-898711-72-1, c LC: QA214.L56. c c Charles Lawson, Richard Hanson, David Kincaid, Fred Krogh, c Basic Linear Algebra Subprograms for FORTRAN usage, c ACM Transactions on Mathematical Software, c Volume 5, Number 3, pages 308-323, 1979. c c Parameters: c c Input, character * ( * ) SRNAME, the name of the routine c which called XERBLA. c c Input, integer INFO, the position of the invalid parameter in the c parameter list of the calling routine. c implicit none integer info character * ( * ) srname write ( *, '(a)' ) ' ' write ( *, '(a)' ) 'XERBLA - Fatal error!' write ( *, '(a)' ) ' On entry to routine ' // trim ( srname ) write ( *, '(a)' ) ' an illegal value was detected for' write ( *, '(a,i6)' ) ' parameter number ', info stop 1 end