1! 2! Copyright (c) 2012-2018, NVIDIA CORPORATION. All rights reserved. 3! 4! Licensed under the Apache License, Version 2.0 (the "License"); 5! you may not use this file except in compliance with the License. 6! You may obtain a copy of the License at 7! 8! http://www.apache.org/licenses/LICENSE-2.0 9! 10! Unless required by applicable law or agreed to in writing, software 11! distributed under the License is distributed on an "AS IS" BASIS, 12! WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. 13! See the License for the specific language governing permissions and 14! limitations under the License. 15! 16 17 18! directives.h -- contains preprocessor directives for F90 rte files 19#include "mmul_dir.h" 20 21subroutine ftn_vmmul_cmplx8( ta, tb, n, k, alpha, a, b, ldb, beta, c ) 22 implicit none 23 integer*8 :: n, k, ldb 24 integer :: ta, tb 25 complex*8, dimension (ldb, * ) :: b 26 complex*8, dimension ( * ) :: a, c 27 complex*8 :: alpha, beta 28 29! local variables 30 integer*8 :: i, j, kk 31 complex*8 :: temp 32 33 34! print *, "#### In vmmul ####" 35 36 if( beta .ne. 0.0 )then 37 do i = 1, n 38 c( i ) = beta * c( i ) 39 enddo 40 else 41 do i = 1, n 42 c( i ) = 0.0 43 enddo 44 end if 45 46 47 48 49 50 if( tb .eq. 2 )then !conjugate b 51 if( ta .eq. 2 )then ! conjugate a - since tb = 2, b is normally oriented 52 if( alpha .eq. ( 1.0, 0.0 ) )then 53 do j = 1, n 54 do kk = 1, k 55 c( j ) = c( j ) + conjg( a( kk ) ) * conjg( b( j, kk ) ) 56 enddo 57 enddo 58 elseif( alpha .eq. (-1.0, 0.0 ) )then 59 do j = 1, n 60 do kk = 1, k 61 c( j ) = c( j ) - conjg( a( kk ) ) * conjg( b( j, kk ) ) 62 enddo 63 enddo 64 else 65 do j = 1, n 66 do kk = 1, k 67 c( j ) = c( j ) + alpha * conjg( a( kk ) ) * conjg( b( kk, j ) ) 68 enddo 69 enddo 70 endif 71 else ! don't conjugate a - if ta != 2, it is just a complex vector 72 if( alpha .eq. ( 1.0, 0.0 ) )then 73 do j = 1, n 74 do kk = 1, k 75 c( j ) = c( j ) + a( kk ) * conjg( b( j, kk ) ) 76 enddo 77 enddo 78 elseif( alpha .eq. ( -1.0, 0.0 ) )then 79 do j = 1, n 80 do kk = 1, k 81 c( j ) = c( j ) - a( kk ) * conjg( b( j, kk ) ) 82 enddo 83 enddo 84 else 85 do j = 1, n 86 do kk = 1, k 87 c( j ) = c( j ) + alpha * a( kk ) * conjg( b( j, kk ) ) 88 enddo 89 enddo 90 endif 91 endif 92 elseif( tb .eq. 1 )then ! b is tranpsosed 93 if( ta .ne. 2 )then ! no conjugation of a is required 94 if( alpha .eq. ( 1.0, 0.0 ) )then 95 do j = 1, n 96 do kk = 1, k 97 c( j ) = c( j ) + a( kk ) * b( j, kk ) 98 enddo 99 enddo 100 elseif( alpha .eq. ( -1.0, 0.0 ) )then 101 do j = 1, n 102 do kk = 1, k 103 c( j ) = c( j ) - a( kk ) * b( j, kk ) 104 enddo 105 enddo 106 else 107 do j = 1, n 108 do kk = 1, k 109 c( j ) = c( j ) + alpha * a( kk ) * b( j, kk ) 110 enddo 111 enddo 112 endif 113 else 114 if( alpha .eq. ( 1.0, 0.0 ) )then 115 do j = 1, n 116 do kk = 1, k 117 c( j ) = c( j ) + conjg( a( kk ) ) * b( j, kk ) 118 enddo 119 enddo 120 elseif( alpha .eq. ( -1.0, 0.0 ) )then 121 do j = 1, n 122 do kk = 1, k 123 c( j ) = c( j ) - conjg( a( kk ) ) * b( j, kk ) 124 enddo 125 enddo 126 else 127 do j = 1, n 128 do kk = 1, k 129 c( j ) = c( j ) + alpha * conjg( a( kk ) ) * b( j, kk ) 130 enddo 131 enddo 132 endif 133 endif 134 else ! b is normally oriented - 135 if( ta .ne. 2 )then ! a is not conjugated 136 if( alpha .eq. ( 1.0, 0.0 ) )then 137 do j = 1, n 138 temp = 0.0 139 do kk = 1, k 140 temp = temp + a(kk) * b( kk, j ) 141 enddo 142 c( j ) = c( j ) + temp 143 enddo 144 elseif( alpha .eq. ( -1.0, 0.0 ) )then 145 do j = 1, n 146 temp = 0.0 147 do kk = 1, k 148 temp = temp + a(kk) * b( kk, j ) 149 enddo 150 c( j ) = c( j ) - temp 151 enddo 152 else 153 do j = 1, n 154 temp = 0.0 155 do kk = 1, k 156 temp = temp + a(kk) * b( kk, j ) 157 enddo 158 c( j ) = c( j ) + alpha * temp 159 enddo 160 endif 161 else ! a is conjugated 162 if( alpha .eq. ( 1.0, 0.0 ) )then 163 do kk = 1, k 164 temp = conjg( a( kk ) ) 165 do j = 1, n 166 c( j ) = c( j ) + temp * b( j, kk ) 167 enddo 168 enddo 169 elseif( alpha .eq. ( -1.0, 0.0 ) )then 170 do kk = 1, k 171 temp = conjg( a( kk ) ) 172 do j = 1, n 173 c( j ) = c( j ) - temp * b( j, kk ) 174 enddo 175 enddo 176 else 177 do kk = 1, k 178 temp = alpha * conjg( a( kk ) ) 179 do j = 1, n 180 c( j ) = c( j ) - temp * b( j, kk ) 181 enddo 182 enddo 183 endif 184 endif 185 endif 186return 187end subroutine ftn_vmmul_cmplx8 188