Actual source code: sgemv_bgl.F

  1: !
  2: !    Fortran kernel for gemv() BLAS operation. This version supports
  3: !  matrix array stored in single precision but vectors in double
  4: !
 5:  #include include/finclude/petscdef.h
  6: !
  7:       subroutine MSGemv_BGL(bs,ncols,A,x,y)
  8:       implicit none
  9:       PetscInt          bs,ncols
 10:       MatScalar        A(bs,ncols)
 11:       PetscScalar      x(ncols),y(bs)

 13:       PetscScalar      bgl_size_0_precision
 14:       PetscInt         i,j

 16:       if( kind(bgl_size_0_precision) .eq. 4) then
 17:         if(sizeof(bgl_size_0_precision) .eq. 4 ) then
 18:            call SGEMV('N', bs, ncols, 1, A, bs, x, 1, 0, y, 1)
 19:            return
 20:         else if(sizeof(bgl_size_0_precision) .eq. 8 ) then
 21:            call CGEMV('N', bs, ncols, 1, A, bs, x, 1, 0, y, 1)
 22:            return
 23:         endif
 24:       else if (kind(bgl_size_0_precision) .eq. 8) then
 25:         if(sizeof(bgl_size_0_precision) .eq. 8 ) then
 26:            call DGEMV('N', bs, ncols, 1, A, bs, x, 1, 0, y, 1)
 27:            return
 28:         else if(sizeof(bgl_size_0_precision) .eq. 16 ) then
 29:            call ZGEMV('N', bs, ncols, 1, A, bs, x, 1, 0, y, 1)
 30:            return
 31:         endif
 32:       endif

 34:       write(6,*)'Error in MSGemv: unavailable type'
 35: 

 37:       do 5, j=1,bs
 38:         y(j) = 0.0d0
 39:  5    continue

 41:       do 10, i=1,ncols
 42:         do 20, j=1,bs
 43:           y(j) = y(j) + A(j,i)*x(i)
 44:  20     continue
 45:  10   continue

 47:       return
 48:       end

 50: 
 51:       subroutine MSGemvp_BGL(bs,ncols,A,x,y)
 52:       implicit none
 53:       PetscInt          bs,ncols
 54:       MatScalar        A(bs,ncols)
 55:       PetscScalar      x(ncols),y(bs)

 57:       PetscInt         i, j
 58:       PetscScalar      bgl_size_0_precision

 60:       if( kind(bgl_size_0_precision) .eq. 4) then
 61:         if(sizeof(bgl_size_0_precision) .eq. 4 ) then
 62:            call SGEMV('N', bs, ncols, 1, A, bs, x, 1, 1, y, 1)
 63:            return
 64:         else if(sizeof(bgl_size_0_precision) .eq. 8 ) then
 65:            call CGEMV('N', bs, ncols, 1, A, bs, x, 1, 1, y, 1)
 66:            return
 67:         endif
 68:       else if (kind(bgl_size_0_precision) .eq. 8) then
 69:         if(sizeof(bgl_size_0_precision) .eq. 8 ) then
 70:            call DGEMV('N', bs, ncols, 1, A, bs, x, 1, 1, y, 1)
 71:            return
 72:         else if(sizeof(bgl_size_0_precision) .eq. 16 ) then
 73:            call ZGEMV('N', bs, ncols, 1, A, bs, x, 1, 1, y, 1)
 74:            return
 75:         endif
 76:       endif
 77: 
 78:       write(6,*)'Error in MSGemvp: unavailable type'

 80:       do 10, i=1,ncols
 81:         do 20, j=1,bs
 82:           y(j) = y(j) + A(j,i)*x(i)
 83:  20     continue
 84:  10   continue


 87:       return
 88:       end


 91:       subroutine MSGemvm_BGL(bs,ncols,A,x,y)
 92:       implicit none
 93:       PetscInt          bs,ncols
 94:       MatScalar        A(bs,ncols)
 95:       PetscScalar      x(ncols),y(bs)

 97:       PetscInt         i, j
 98:       PetscScalar      bgl_size_0_precision

100:       if( kind(bgl_size_0_precision) .eq. 4) then
101:         if(sizeof(bgl_size_0_precision) .eq. 4 ) then
102:            call SGEMV('N', bs, ncols, 1, A, bs, x, 1, -1, y, 1)
103:            return
104:         else if(sizeof(bgl_size_0_precision) .eq. 8 ) then
105:            call CGEMV('N', bs, ncols, 1, A, bs, x, 1, -1, y, 1)
106:            return
107:         endif
108:       else if (kind(bgl_size_0_precision) .eq. 8) then
109:         if(sizeof(bgl_size_0_precision) .eq. 8 ) then
110:            call DGEMV('N', bs, ncols, 1, A, bs, x, 1, -1, y, 1)
111:            return
112:         else if(sizeof(bgl_size_0_precision) .eq. 16 ) then
113:            call ZGEMV('N', bs, ncols, 1, A, bs, x, 1, -1, y, 1)
114:            return
115:         endif
116:       endif

118:       write(6,*)'Error in MSGemvm: unavailable type'

120:       do 10, i=1,ncols
121:         do 20, j=1,bs
122:           y(j) = y(j) - A(j,i)*x(i)
123:  20     continue
124:  10   continue

126:       return
127:       end


130:       subroutine MSGemvt_BGL(bs,ncols,A,x,y)
131:       implicit none
132:       PetscInt          bs,ncols
133:       MatScalar        A(bs,ncols)
134:       PetscScalar      x(bs),y(ncols)

136:       PetscInt          i,j
137:       PetscScalar      sum
138:       PetscScalar      bgl_size_0_precision

140:       if( kind(bgl_size_0_precision) .eq. 4) then
141:         if(sizeof(bgl_size_0_precision) .eq. 4 ) then
142:            call SGEMV('T', bs, ncols, 1, A, bs, x, 1, y, 1)
143:            return
144:         else if(sizeof(bgl_size_0_precision) .eq. 8 ) then
145:            call CGEMV('T', bs, ncols, 1, A, bs, x, 1, y, 1)
146:            return
147:         endif
148:       else if (kind(bgl_size_0_precision) .eq. 8) then
149:         if(sizeof(bgl_size_0_precision) .eq. 8 ) then
150:            call DGEMV('T', bs, ncols, 1, A, bs, x, 1, y, 1)
151:            return
152:         else if(sizeof(bgl_size_0_precision) .eq. 16 ) then
153:            call ZGEMV('T', bs, ncols, 1, A, bs, x, 1, y, 1)
154:            return
155:         endif
156:       endif

158:       write(6,*)'Error in MSGemvt: unavailable type'
159: 
160:       do 10, i=1,ncols
161:         sum = y(i)
162:         do 20, j=1,bs
163:           sum = sum + A(j,i)*x(j)
164:  20     continue
165:         y(i) = sum
166:  10   continue

168:       return
169:       end

171:       subroutine MSGemm_BGL(bs,A,B,C)
172:       implicit none
173:       PetscInt    bs
174:       MatScalar   A(bs,bs),B(bs,bs),C(bs,bs)
175:       PetscScalar sum
176:       PetscInt    i,j,k

178:       PetscScalar bgl_size_0_precision

180:       if( kind(bgl_size_0_precision) .eq. 4) then
181:         if(sizeof(bgl_size_0_precision) .eq. 4 ) then
182:            call SGEMM('N', 'N', bs, bs, bs, -1, B, bs, C, bs, 1, A, bs)
183:            return
184:         else if(sizeof(bgl_size_0_precision) .eq. 8 ) then
185:            call CGEMM('N', 'N', bs, bs, bs, -1, B, bs, C, bs, 1, A, bs)
186:            return
187:         endif
188:       else if (kind(bgl_size_0_precision) .eq. 8) then
189:         if(sizeof(bgl_size_0_precision) .eq. 8 ) then
190:            call DGEMM('N', 'N', bs, bs, bs, -1, B, bs, C, bs, 1, A, bs)
191:            return
192:         else if(sizeof(bgl_size_0_precision) .eq. 16 ) then
193:            call ZGEMM('N', 'N', bs, bs, bs, -1, B, bs, C, bs, 1, A, bs)
194:            return
195:         endif
196:       endif


199:       write(6,*)'Error in MSGemm: unavailable type'

201:       do 10, i=1,bs
202:         do 20, j=1,bs
203:           sum = A(i,j)
204:           do 30, k=1,bs
205:             sum = sum - B(i,k)*C(k,j)
206:  30       continue
207:           A(i,j) = sum
208:  20     continue
209:  10   continue

211:       return
212:       end


215:       subroutine MSGemmi_BGL(bs,A,C,B)
216:       implicit none
217:       PetscInt    bs
218:       MatScalar   A(bs,bs),B(bs,bs),C(bs,bs)
219:       PetscScalar sum

221:       PetscInt    i,j,k

223:       PetscScalar bgl_size_0_precision

225:       if( kind(bgl_size_0_precision) .eq. 4) then
226:         if(sizeof(bgl_size_0_precision) .eq. 4 ) then
227:            call SGEMM('N', 'N', bs, bs, bs, 1, B, bs, C, bs, 0, A, bs)
228:            return
229:         else if(sizeof(bgl_size_0_precision) .eq. 8 ) then
230:            call CGEMM('N', 'N', bs, bs, bs, 1, B, bs, C, bs, 0, A, bs)
231:            return
232:         endif
233:       else if (kind(bgl_size_0_precision) .eq. 8) then
234:         if(sizeof(bgl_size_0_precision) .eq. 8 ) then
235:            call DGEMM('N', 'N', bs, bs, bs, 1, B, bs, C, bs, 0, A, bs)
236:            return
237:         else if(sizeof(bgl_size_0_precision) .eq. 16 ) then
238:            call ZGEMM('N', 'N', bs, bs, bs, 1, B, bs, C, bs, 0, A, bs)
239:            return
240:         endif
241:       endif

243:       write(6,*)'Error in MSGemmi: unavailable type'

245:       do 10, i=1,bs
246:         do 20, j=1,bs
247:           sum = 0.0d0
248:           do 30, k=1,bs
249:             sum = sum + B(i,k)*C(k,j)
250:  30       continue
251:           A(i,j) = sum
252:  20     continue
253:  10   continue

255:       return
256:       end