c********************************************************** c * c TEST SUITE FOR VECTORIZING COMPILERS * c * c Version: 2.0 * c Date: 3/14/88 * c Authors: Original loops from a variety of * c sources. Collection/synthesis by: * c * c David Callahan - Rice University * c Jack Dongarra - Argonne National Lab * c David Levine - Argonne National Lab * c * c Compile: Use the compiler options that show the * c automatic vectorization capabilities * c of the compiler. Do not use compiler * c directives or interactive compilation * c features for additional vectorizations. * c Do not modify the source code. The * c subroutines in section 1.5 may be * c compiled independently. * c * c********************************************************** c********************************************************** c * c DEPENDENCE ANALYSIS * c * c********************************************************** c%1.1 subroutine s111(a,n) c c linear dependence testing c dependence testing - vectorizable c dimension a(1000) do 400 i = 2,100,2 a(i) = a(i-1) 400 continue write (unit=6,fmt=100) a(n) 100 format (e12.6) return end c%1.1 subroutine s112(a,n) c c linear dependence testing c loop reversal c dimension a(1000) do 700 i = 999,1,-1 a(i+1) = a(i) 700 continue write (unit=6,fmt=100) a(n) 100 format (e12.6) return end c%1.1 subroutine s113(a,n) c c linear dependence testing c a(i)=a(1) but no actual dependence cycle c integer n real a(*) do 610 i = 2,n a(i) = a(1) 610 continue return end c%1.1 subroutine s114(aa,bb,n) c c linear dependence testing c transpose vectorization c integer n real aa(n,*),bb(n,*) do 300 j = 1,n do 300 i = 1,j-1 aa(i,j) = aa(j,i) + bb(i,j) 300 continue return end c%1.1 subroutine s115(aa,n) c c linear dependence testing c lower triangular system c c integer n real aa(n,*) do 320 j = 1,n do 320 k = 1,j-1 do 320 i = k+1,n 320 aa(i,j) = aa(i,j) + aa(i,k) * aa(k,j) return end c%1.1 subroutine s116(a,b,n) c c linear dependence testing c integer n real a(*),b(*) do 450 i = 1,n,5 a(i) = a(i+1) + a(i)*b(i) a(i+1) = a(i+2) + a(i+1)*b(i+1) a(i+2) = a(i+3) + a(i+2)*b(i+2) a(i+3) = a(i+4) + a(i+3)*b(i+3) a(i+4) = a(i+5) + a(i+4)*b(i+4) 450 continue return end c%1.1 subroutine s117(a,logn,n) c c linear dependence testing c fft subscripting c integer n,logn real a(*) j = 1 do 970 i = 1,logn do 971 k = 1,n,j do 971 l = 1,j t = a(k+l-1) + a(k+l+j-1) u = a(k+l-1) - a(k+l+j-1) a(k+l-1) = t a(k+l+j-1) = u 971 continue j = j * 2 970 continue return end c%1.2 subroutine s121 c c induction variable recognition c loop with possible ambiguity because of scalar store c parameter ( n = 1000 ) real a(1000) integer i, j do 70, i = 1, n-1 j = i+1 a(i) = a(j) 70 continue write(unit=6,fmt=11) (a(i),i=1,n) 11 format(e12.6) return end c%1.2 subroutine s122(xx,yy,n) c c induction variable recognition c mixed variable and constant lb, ub, and stride c@ k is not initialized before first use. c@ n1,n3 are not defined. c dimension xx(100), yy(100), zz(100) xx(1) = 6. j = 6 do 60 i=n1,100,n3 xx(i)=yy(i)*zz(100-k+1) k = k + j 60 continue write (unit=6,fmt=11) j, xx(n) 11 format (i6,e12.6) return end c%1.2 subroutine s123(a,b,c,n) c c induction variable under an if c integer n real a(*),b(*),c(*) j = 0 do 50 i = 1,n j = j + 1 a(j) = b(i) if(c(i).gt.0) then j = j + 1 a(j) = c(i) endif 50 continue return end c%1.2 subroutine s124(a,b,c,n) c c induction variable recognition c induction variable under both sides of if (same value) c integer n real a(*),b(*),c(*) j = 0 do 60 i = 1,n if(b(i).gt.0) then j = j + 1 a(j) = b(i) else j = j + 1 a(j) = c(i) endif 60 continue return end c%1.2 subroutine s125(a,bb,cc,n) c c induction variable recognition c induction variable in two loops; collapsing possible c integer n real a(*),bb(n,*),cc(n,*) k = 0 do 20 i = 1,n do 20 j = 1,n k = k + 1 a(k) = bb(i,j) + cc(i,j) 20 continue return end c%1.2 subroutine s126(a,bb,n) c c induction variable recognition c induction variable in two loops; recurrence in inner loop c integer n real a(*),bb(n,*) k = 1 do 40 i = 1,n do 41 j = 2,n bb(i,j) = bb(i,j-1) + a(k) k = k + 1 41 continue k = k + 1 40 continue return end c%1.2 subroutine s127(a,b,c,n) c c induction variable recognition c induction variable with multiple increments c integer n real a(*),b(*),c(*) j = 0 do 10 i = 1,n j = j + 1 a(j) = b(i) j = j + 1 a(j) = c(i) 10 continue return end c%1.3 subroutine s131(a,b,n) c c global data flow analysis c forward substitution c integer n real a(*),b(*) m = 1 if(a(1).gt.0)then a(1) = b(1) endif do 340 i = 2,n-1 a(i) = a(i+m) + b(i) 340 continue return end c%1.3 subroutine s132 c c global data flow analysis c loop with multiple dimension ambiguous subscripts c parameter ( n = 100 ) real a(100,100) integer i, j, k, m m = 1 j = m k = m+1 do 170, i=2,n a(i,j) = a(i-1,k)*3.5 170 continue write(unit=6,fmt=11) ((a(i,j),i=1,n),j=1,n) 11 format(e12.6) return end c%1.5 subroutine s151(a,b,n) c c interprocedural data flow analysis c universal compilation - passing parameter c information into a subroutine c integer n real a(*),b(*) call s151s(a,b,n-1,1) return end subroutine s151s(a,b,n,m) integer n,m real a(*),b(*) do 730 i = 1,n a(i) = a(i+m) + b(i) 730 continue return end c%1.5 subroutine s152(a,b,n) c c interprocedural data flow analysis c universal compilation - collecting info from subroutine c integer n real a(*),b(*) do 750 i = 1,n b(i) = b(i) + 2 call s152s(a,b,i) 750 continue return end subroutine s152s(a,b,i) integer i real a(*),b(*) a(i) = a(i) + b(i) return end c%1.6 subroutine s161(a,b,x,n) c c control flow c tests for recognition of loop independent dependences c between statements in mutually exclusive regions. c@ array bounds error when i = 100. c dimension a(100),b(100),x(100) do 30 i = 1,100 if (b(i).lt.0) go to 10 a(i) = x(i) go to 30 10 x(i+1) = a(i) 30 continue write(unit=6,fmt=11) a(n),b(n),x(n) 11 format(3e12.6) return end c%1.7 subroutine s171(a,b,n) c c symbolics c symbolic dependence tests c integer n real a(*),b(*) do 1030 i = 1,n a(i*n) = a(i*n) + b(i) 1030 continue return end c%1.7 subroutine s172(a,b,n) c c symbolics c symbolics - not vectorizable c@ could be vectorized if you assume n3 .NE. 0, c@ a reasonable assumption. c@ n1,n2,n3 are not defined. c dimension a(1000),b(1000) do 500 i = n1,n2,n3 a(i) = a(i) + b(i) 500 continue write (unit=6,fmt=100) a(n) 100 format (e12.6) return end c%1.7 subroutine s173(a,b,n) c c symbolics c expression in loop bounds and subscripts c integer n real a(*),b(*) do 370 i = 1,n/2 a(i+n/2) = a(i) + b(i) 370 continue return end c%1.7 subroutine s174 c c symbolics c loop with subscript that may seem ambiguous c parameter ( n = 1000 ) real a(1000) integer i do 50, i= 1, n/2 a(i) = a( i + n/2 ) 50 continue write(unit=6,fmt=11) (a(i),i=1,n) 11 format(e12.6) return end c%1.7 subroutine s175(a,b,n,inc) c c symbolics c symbolic dependence tests c integer n,inc real a(*),b(*) do 1020 i = 1,n,inc a(i) = a(i+inc) + b(i) 1020 continue return end c c********************************************************** c * c VECTORIZATION * c * c********************************************************** c%2.1 subroutine s211(a,b,c,n) c c statement reordering c statement reordering allows vectorization c integer n real a(*),b(*),c(*) do 270 i = 2,n-1 a(i) = b(i-1) + c(i) b(i) = b(i+1) - 2. 270 continue return end c%2.1 subroutine s212 c c statement reordering c dependency needing temporary c parameter ( n = 1000 ) real a(1000), b(1000), x(6000) integer i do 20, i=1,n-1 a(i) = x(i) b(i) = b(i)+a(i+1) 20 continue write(unit=6,fmt=11) (a(i),i=1,n) write(unit=6,fmt=11) (b(i),i=1,n) 11 format(e12.6) return end c%2.2 subroutine s221 c c loop distribution c loop that is partially recursive c parameter ( n = 1000 ) real a(1000), b(1000), x(6000), y(6000) integer i do 80, i = 2 , n a(i) = a(i) + ( x(i) * y(i) ) b(i) = b(i-1) + a(i) + y(i) 80 continue write(unit=6,fmt=11) (a(i),i=1,n) write(unit=6,fmt=11) (b(i),i=1,n) 11 format(e12.6) return end c%2.2 subroutine s222(a,b,c,n) c c loop distribution c partial loop vectorization, recurrence in middle c integer n real a(*),b(*),c(*) do 240 i = 2,n a(i) = a(i) + b(i) b(i) = b(i-1)*b(i-1)*a(i) a(i) = a(i) - b(i) 240 continue return end c%2.3 subroutine s231 c c loop interchange c loop with multiple dimension recursion c parameter ( n = 100 ) real a(100,100), b(100,100) integer i, j do 160, i=1,n do 160, j=2,n a(i,j) = a(i,j-1)+b(i,j) 160 continue write(unit=6,fmt=11) ((a(i,j),i=1,n),j=1,n) 11 format(e12.6) return end c%2.3 subroutine s232(aa,bb,n) c c loop interchange c interchanging of triangular loops c integer n real aa(n,*),bb(n,*) do 290 j = 2,n do 290 i = 2,j aa(i,j) = aa(i-1,j)*aa(i-1,j)+bb(i,j) 290 continue return end c%2.3 subroutine s233(aa,bb,cc,n) c c loop interchange c interchanging with one of two inner loops c integer n real aa(n,*),bb(n,*),cc(n,*) do 840 i = 2,n do 841 j = 2,n aa(i,j) = aa(i,j-1) + cc(i,j)/aa(i,j-1) 841 continue do 842 j = 2,n bb(i,j) = bb(i-1,j) + cc(i,j)/bb(i-1,j) 842 continue 840 continue return end c%2.3 subroutine s234(aa,bb,cc,n) c c loop interchange c if loop to do loop, interchanging with if loop necessary c integer n real aa(n,*),bb(n,*),cc(n,*) i = 1 232 if(i.gt.n) goto 231 j = 3 233 if(j.gt.n) goto 230 aa(i,j) = aa(i,j-1)*bb(i,j-1) bb(i,j) = aa(i,j-1)+bb(i,j-2) cc(i,j) = cc(i,j-1)*cc(i,j-2) j = j + 1 goto 233 230 i = i + 1 goto 232 231 continue return end c%2.4 subroutine s241(a,b,c,n) c c node splitting c preloading necessary to allow vectorization c integer n real a(*),b(*),c(*) do 280 i = 2,n-1 a(i) = b(i) + c(i) b(i) = a(i) + a(i+1) 280 continue return end c%2.4 subroutine s242(a,b,c,n) c c node splitting c@ array bounds error when i=1, i=100 c dimension a(100),b(100),c(100),d(100) dimension abc1(100),abc2(100),bcd1(100) do 510 i = 1,100 abc1(i-1) = abc2(i-1) + 1 bcd1(i-1) = abc1(i-1) ** 2 a(i) = bcd1(i-1) + d(i) abc2(i) = b(i+1) + b(i-1) b(i) = c(i) + 1 c(i+1) = b(i) + 1 510 continue write(unit=6,fmt=12) c(n) write(unit=6,fmt=12) b(n) 12 format (e12.6) return end c%2.4 subroutine s243(a,b,c,n) c c node splitting c false dependence cycle breaking c integer n real a(*),b(*),c(*) do 630 i = 2,n-1 a(i) = b(i) + c(i) b(i) = a(i) - 2. a(i) = b(i) + a(i+1) 630 continue return end c%2.4 subroutine s244(a,b,c,n) c c node splitting c false dependence cycle breaking c integer n real a(*),b(*),c(*) do 640 i = 2,n-1 a(i) = b(i) + c(i) b(i) = c(i) - 2. a(i+1) = b(i) + a(i+1) 640 continue return end c%2.4 subroutine s245(a,b,c,n) c c node splitting c array bounds error when i=1 c dimension a(100),b(100),c(100),d(100) t2 = 0 do 800 i = 1,100 a(i) = a(i-1) + t1 + t2 + b(i) + c(i) + d(i) 800 continue write(unit=6,fmt=12) a(n) 12 format (e12.6) return end c%2.5 subroutine s251(b,c,n) c c scalar and array expansion c scalar expansion c dimension b(1000),c(1000),d(1000) do 900 i = 1,100,1 abc = b(i) * c(i) d(i) = abc * abc 900 continue write (unit=6,fmt=100) d(n) 100 format (e12.6) return end c%2.5 subroutine s252 c c scalar and array expansion c loop with ambiguous scalar temporary c parameter ( n = 6000 ) real a(6000), x(6000), y(6000), s, t integer i t = 0. do 40, i=1,n s = x(i) * y(i) a(i) = s + t t = s 40 continue write(unit=6,fmt=11) (a(i),i=1,n) 11 format(e12.6) return end c%2.5 subroutine s253(a,b,c,n) c c scalar and array expansion c scalar expansion, assigned under if c integer n real a(*),b(*),c(*) do 760 i = 1,n if(a(i).gt.b(i))then t = a(i) - b(i) c(i) = c(i) + t a(i) = t endif 760 continue return end c%2.5 subroutine s254(a,b,n) c c scalar and array expansion c carry around variable c integer n real a(*),b(*) x = a(n) do 390 i = 1,n b(i) = (a(i) + x) / 2. x = a(i) 390 continue return end c%2.5 subroutine s255(a,b,n) c c scalar and array expansion c carry around variables, 2 levels c integer n real a(*),b(*) x = a(n) y = a(n-1) do 400 i = 1,n b(i) = (a(i) + x + y) / 3. y = x x = a(i) 400 continue return end c%2.7 subroutine s271 c c control flow c loop with singularity handling c parameter ( n = 1000 ) real a(1000), y(6000), z(6000) integer i do 140, i=1,n if (y(i).gt.0.) a(i) = y(i)/z(i) 140 continue write(unit=6,fmt=11) (a(i),i=1,n) 11 format(e12.6) return end c%2.7 subroutine s272 c c control flow c loop with independent conditional c parameter ( n = 1000 ) real a(1000), x(6000), y(6000), z(6000) real s, r, t integer i t = 1. do 100, i = 1, n if (z(i) .ge. t) then s = x(i) * y(i) + 3.1 r = x(i) + y(i) * 2.9 a(i) = sqrt(s**2*r) endif 100 continue write(unit=6,fmt=11) (a(i),i=1,n) 11 format(e12.6) return end c%2.7 subroutine s273 c c control flow c simple loop with dependent conditional c parameter ( n = 4000 ) real a(4000), b(4000), c(4000) real x(6000), y(6000), z(6000) integer i do 120, i = 1 , n a(i) = a(i) + y(i) + z(i) if (a(i) .lt. 0.) b(i) = b(i) + x(i) + y(i) c(i) = c(i) + a(i) + x(i) 120 continue write(unit=6,fmt=11) (a(i),i=1,n) write(unit=6,fmt=11) (b(i),i=1,n) write(unit=6,fmt=11) (c(i),i=1,n) 11 format(e12.6) return end c%2.7 subroutine s274 c c control flow c complex loop with dependent conditional c parameter ( n = 1000 ) real a(1000), b(1000) real x(6000), y(6000), z(6000) integer i do 130, i = 1 , n a(i) = x(i) + z(i) if (a(i) .eq. 0.) then b(i) = a(i) * b(i) else a(i) = y(i) * z(i) b(i) = 1. endif 130 continue write(unit=6,fmt=11) (a(i),i=1,n) write(unit=6,fmt=11) (b(i),i=1,n) 11 format(e12.6) return end c%2.7 subroutine s275(aa,bb,n) c c control flow c if around inner loop, interchanged needed c integer n real aa(n,*),bb(n,*) do 480 i = 2,n if(bb(i,1).gt.0)then do 481 j = 2,n 481 bb(i,j) = bb(i,j)/bb(i,j-1) endif 480 continue return end c%2.7 subroutine s276(a,b,c,n) c c control flow c ifs which trim the index set c integer n real a(*),b(*),c(*) do 660 i = 1,n a(i) = b(i) + c(i) if(i.gt.5) b(i) = abs(b(i)) if(i.lt.99) a(i) = -a(i) 660 continue return end c%2.7 subroutine s277(a,b,x,n) c c control flow c test for dependences arising from guard variable computation. c@ y is used before being defined c@ array bounds error when i = 100. c dimension a(100),b(100),x(100),y(100) do 130 i = 1,100 if (a(i).ge.0) go to 120 if (b(i).ge.0) go to 110 a(i) = x(i) 110 continue b(i+1) = y(i) 120 continue 130 continue write(unit=6,fmt=11) a(n),b(n),x(n) 11 format(3e12.6) return end c%2.7 subroutine s278(a,b,c,n) c c control flow c if/goto to block if-then-else c integer n real a(*),b(*),c(*) do 120 i = 1,n if(a(i).gt.0)goto 121 b(i) = -b(i) goto 122 121 continue c(i) = -c(i) 122 continue a(i) = b(i) + c(i) 120 continue return end c%2.7 subroutine s279(a,b,c,n) c c control flow c vector if/gotos c integer n real a(*),b(*),c(*) do 810 i = 1,n if(a(i).gt.0)goto 811 b(i) = -b(i) if(abs(b(i)).le.a(i))goto 812 c(i) = abs(c(i)) goto 812 811 continue c(i) = -c(i) 812 continue a(i) = b(i) + c(i) 810 continue return end c%2.7 subroutine s2710(a,b,c,x,n) c c control flow c scalar and vector ifs c integer n real a(*),b(*),c(*),x do 790 i = 1,n if(a(i).gt.b(i))then a(i) = a(i) - b(i) if(n.gt.10)then c(i) = abs(c(i)) else c(i) = 0. endif else b(i) = a(i) if(x.gt.0)then c(i) = a(i) else c(i) = -c(i) endif endif 790 continue return end c%2.7 subroutine s2711(a,b,c,n) c c control flow c semantic if removal c integer n real a(*),b(*),c(*) do 650 i = 1,n if(a(i).ne.0) b(i) = b(i) + a(i) * c(i) 650 continue return end c%4.5 subroutine s2712(a,b,n) c c control flow c if to elemental min c integer n real a(*),b(*) do 70 i = 1,n if(a(i).gt.b(i)) a(i) = b(i) 70 continue return end c%2.8 subroutine s281(a,b,n) c c crossing thresholds ( index set splitting ) c index set splitting c integer n real a(*),b(*) do 990 i = 1,n x = a(n-i+1) + b(i) a(i) = x + 3. b(i) = x 990 continue return end c%2.9 subroutine s291(a,b,n) c c loop peeling c wrap around variable, 1 level c integer n real a(*),b(*) im1 = n do 410 i = 1,n b(i) = (a(i) + a(im1)) / 2. im1 = i 410 continue return end c%2.9 subroutine s292(a,b,n) c c loop peeling c wrap around variable, 2 levels c integer n real a(*),b(*) im1 = n im2 = n-1 do 420 i = 1,n b(i) = (a(i) + a(im1) + a(im2)) / 3. im2 = im1 im1 = i 420 continue return end c%2.9 subroutine s293(a,n) c c loop peeling c a(i)=a(1) with actual dependence cycle c@ above comment misleading, loop is vectorizable c integer n real a(*) do 620 i = 1,n a(i) = a(1) 620 continue return end c c********************************************************** c * c IDIOM RECOGNITION * c * c********************************************************** c%3.1 subroutine s311(a,b,x,n) c c reductions c sum reduction c integer n real a(*),b(*),x do 850 i = 1,n x = x + a(i) b(i) = a(i) + 2. 850 continue return end c%3.1 subroutine s312(a,b,x,n) c c reductions c product reduction c integer n real a(*),b(*),x do 860 i = 1,n x = x * a(i) b(i) = a(i) + 2. 860 continue return end c%3.1 subroutine s313(a,b,n) c c reductions c dot product c dimension a(1000),b(1000) s = 0. do 930 i = 1,n s = s + a(i) * b(i) 930 continue write (unit=6,fmt=100) s 100 format (e12.6) return end c%3.1 subroutine s314(x,b,n) c c reductions c if to max reduction c integer n real x,b(*) do 80 i = 1,n if(b(i).gt.x) x = b(i) 80 continue return end c%3.1 subroutine s315(x,j,b,n) c c reductions c if to max with index reduction, 1 dimension c integer n,j real x,b(*) do 90 i = 1,n if(b(i).gt.x)then x = b(i) j = i endif 90 continue return end c%3.1 subroutine s316(a,n) c c reductions c minval c dimension a(1000) s = a(1) do 960 i = 2,n if (a(i) .lt. s) s = a(i) 960 continue write (unit=6,fmt=100) s write (unit=6,fmt=100) a(n) 100 format (e12.6) return end c%3.1 subroutine s317(n) c c reductions c tests scalar expansion. From a benchmark to test the c scalar speed of a machine. The best results are fully c vectorized. the compiler expands .995 into a vector of c .995 and does a product reduction on the vector of .995. c@ note: this loop has closed form solution: q = .995**n c q = 1. do 50, i = 1,n q = .995*q 50 continue write(unit=6,fmt=11) q 11 format(e12.6) return end c%3.1 subroutine s318(a,b,x,n) c c reductions c integer n real a(*),b(*),x do 430 i = 1,n,5 x = x + a(i)*b(i) + a(i+1)*b(i+1) + a(i+2)*b(i+2) $ + a(i+3)*b(i+3) + a(i+4)*b(i+4) 430 continue return end c%3.2 subroutine s321(a,b,n) c c recurrences c first order linear recurrence c integer n real a(*),b(*) do 870 i = 2,n a(i) = a(i) + a(i-1)*b(i) 870 continue return end c%3.2 subroutine s322(a,b,c,n) c c recurrences c second order linear recurrence c integer n real a(*),b(*),c(*) do 880 i = 3,n a(i) = a(i) + a(i-1)*b(i) + a(i-2)*c(i) 880 continue return end c%3.2 subroutine s323(a,b,c,d,n) c c recurrences c coupled recurrence c@ array bounds error when i = 1. c integer n real a(*),b(*),c(*),d(*) do 1040 i = 1,n a(i) = b(i-1) + c(i) b(i) = a(i) + d(i) 1040 continue return end c%3.3 subroutine s331(a,j,n) c c search loops c if to last-1 c integer j,n real a(*) do 130 i = 1,n if(a(i).eq.0) j = i 130 continue return end c%3.3 subroutine s332(a,b,j,n) c c search loops c search loop saving index c integer n,j real a(*),b(*) do 510 i = 1,n if(a(i).gt.b(i))then j = i goto 511 endif 510 continue j = 0 return 511 continue return end c%3.4 subroutine s341(a,b,j,n) c c packing c integer n,j real a(*),b(*) j = 0 do 900 i = 1,n if(a(i).gt.0)then j = j + 1 b(j) = a(i) endif 900 continue return end c%3.4 subroutine s342(a,b,j,n) c c packing c unpacking c integer n,j real a(*),b(*) j = 0 do 910 i = 1,n if(a(i).gt.0)then j = j + 1 a(i) = b(j) endif 910 continue return end c c********************************************************** c * c LANGUAGE COMPLETENESS * c * c********************************************************** c%4.1 subroutine s411(a,b,c,n) c c loop recognition c if loop to do loop, zero trip c integer n real a(*),b(*),c(*) i = 0 140 continue i = i + 1 if(i.gt.n)goto 141 a(i) = b(i) + c(i) goto 140 141 continue return end c%4.1 subroutine s412(a,b,c,n,inc) c c loop recognition c if loop with variable increment c integer n,inc real a(*),b(*),c(*) i = 0 950 continue i = i + inc if(i.gt.n)goto 951 a(i) = b(i) + c(i) goto 950 951 continue return end c%4.1 subroutine s413(a,b,c,n) c c loop recognition c if loop to do loop, code on both sides of increment c integer n real a(*),b(*),c(*) i = 0 180 continue if(i.gt.n)goto 181 b(i) = abs(b(i)) i = i + 1 a(i) = c(i) goto 180 181 continue return end c%4.1 subroutine s414(aa,bb,cc,n) c c loop recognition c if loop to do loop, interchanging with do necessary c integer n real aa(n,*),bb(n,*),cc(n,*) i = 1 222 if(i.gt.n) goto 221 do 220 j = 3,n aa(i,j) = aa(i,j-1)*bb(i,j-1) bb(i,j) = aa(i,j-1)+bb(i,j-2) 220 cc(i,j) = cc(i,j-1)*cc(i,j-2) i = i + 1 goto 222 221 continue return end c%4.2 subroutine s421(n) c c storage classes and equivalencing c equivalence- no overlap c@ array bounds error when i = 100. c dimension xx(100), yy(100) equivalence (xx(1),yy(1)) do 205 i = 1,100 xx(i) = yy(i+1) 205 continue write(unit=6,fmt=100) xx(n) 100 format(e12.6) return end c%4.2 subroutine s422(n) c c storage classes and equivalencing c equivalence- vectorizable c@ array bounds error when i,j c@ are simultaneously > 96. c dimension x(100,100),y(100,100) equivalence (x(1,1),y(1,1)) do 410 i = 1,100 do 420 j = 1,100 x(j+4,i+4) = y(j,i) 420 continue 410 continue write(unit=6,fmt=100) x(n,n) 100 format(e12.6) return end c%4.2 subroutine s423(n) c c storage classes and equivalencing c common and equivalenced variables - no overlap c@ misleading comment, there is an anti-dependence c dimension cc(100) common /com1/ii1 common /com2/aa(200) equivalence (aa(50),cc(1)) do 116 ii1 = 1,100 aa(ii1+1) = cc(ii1) 116 continue write(unit=6,fmt=100) aa(n) 100 format(e12.6) return end c%4.2 subroutine s424(n) c c storage classes and equivalencing c common and equivalenced variables - no overlap c@ threshold is 100 >= loop upper bound => no dependence c dimension cc(100) common /com1/ii1 common /com2/aa(200) equivalence (aa(50),cc(1)) do 118 ii1 = 1,100 cc(ii1+50) = aa(ii1-1) 118 continue write(unit=6,fmt=100) cc(n) 100 format(e12.6) return end c%4.2 subroutine s425(n) c c storage classes and equivalencing c common and equivalence statement c@ anti-dependence with threshold of 4 c dimension b(1000) c common /comm1/cc1(100),cc2(100) common /comm1/cc1(100),cc2(100),cc3(100),cc4(100,100) dimension eqv1(100),eqv2(90) equivalence (eqv1(1),cc2(1)) equivalence (eqv2(1),cc1(5)) do 920 i = 1,85 eqv2(i) = cc1(i+8) + b(i) 920 continue write(6,100) eqv2(n) 100 format(e12.6) return end c%4.2 subroutine s426(n,y) c c storage classes and equivalencing c common and equivalence statement c@ anti dependence with distance vector <2,6> c@ vectorizable with respect to both loops c dimension y(100,100) common /comm1/cc1(100),cc2(100),cc3(100),cc4(100,100) dimension eqv1(100),eqv2(90),eqv3(100),eqv4(100,99) equivalence (eqv1(1),cc2(1)) equivalence (eqv2(1),cc1(5)) equivalence (eqv3(5),cc3(5)) equivalence (eqv4(1,2),cc4(1,1)) do 940 j = 2,80 do 950 i = 1,90 eqv4(i,j) = cc4(i+2,j+5) + y(i,j) 950 continue 940 continue write(6,100) eqv4(n,n) 100 format(e12.6) return end c%4.2 subroutine s427(n) c c storage classes and equivalencing c common and equivalenced variables - overlap c@ a partially negative test, vectorizable in chunks of <=50 c dimension c(100) common /com1/i common /com2/a(200) equivalence (a(50),c(1)) do 115 i = 1,100 c(i+1) = a(i) 115 continue write(unit=6,fmt=100) c(n) 100 format(e12.6) return end c%4.4 subroutine s441(a,c,n) c c non-logical if's c arithmetic if c@ xx is used before being defined. c dimension a(1000),b(1000),c(1000),d(1000),xx(100) do 710 i = 1,100 if (d(i)) 720,730,740 720 c(i) = a(i) goto 750 730 c(i) = b(i) goto 750 740 c(i) = xx(i) 750 continue 710 continue write(6,100) c(n) 100 format(e12.6) return end c%4.4 subroutine s442(a,c,n) c c non-logical if's c computed goto c@ aa is used before being defined. c@ ii is used before being defined. c dimension a(1000),b(1000),c(1000),d(1000),aa(200),ii(1000) do 810 i = 1,100 goto (815,820,830,840) ii(i) 815 c(i) = aa(i) goto 850 820 c(i) = a(i) goto 850 830 c(i) = b(i) goto 850 840 c(i) = d(i) 850 continue 810 continue write(6,100) c(n) 100 format(e12.6) return end c%4.5 subroutine s451(a,b,c,n) c c intrinsic functions c intrinsics c integer n real a(*),b(*),c(*) do 930 i = 1,n a(i) = sin(b(i)) + cos(c(i)) 930 continue return end c%4.5 subroutine s452(a,b,c,m) c c intrinsic functions c seq function c dimension a(1000),b(1000),c(1000) do 1250 i = 1,m a(i) = b(i) + c(i) + i 1250 continue return end c%4.6 subroutine s461(a,c,n) c c i/o statements c write statement c dimension a(1000),b(1000),c(1000),d(1000) do 600 i = 1,100 b(i) = d(i) c(i) = a(i) + b(i) write(6,100) c(i) 600 continue write(6,100) b(n) 100 format(e12.6) return end c%4.7 subroutine s471(z,n) c c call statements c dimension x(100,100),z(100,100),w(100,100),y(100,100) equivalence (x(1,1),y(1,1)) do 510 i = 1,100 do 520 j = 1,100 x(j,i) = z(j,i) call sub2 z(j,i) = w(j,i) 520 continue 510 continue write(unit=6,fmt=100) z(n,n) 100 format(e12.6) return end c%4.8 subroutine s481(a,b,c,n) c c non-local goto's c stop statement c dimension a(1000),b(1000),c(1000) do 110 i = 1,100 if (a(i) .lt. 0.) stop 'stop 1' b(i) = c(i) 110 continue write (unit=6,fmt=100) b(n) 100 format(e12.6) return end c%4.8 subroutine s482(a,b,c,n) c c non-local goto's c other loop exit with code before exit c integer n real a(*),b(*),c(*) do 520 i = 1,n c(i) = a(i) if(a(i).gt.b(i))goto 521 520 continue return 521 continue return end c%4.10 subroutine s4101(a,b,c,n) c c data types c complex arithmetic c integer n complex a(*),b(*),c(*) do 920 i = 1,n c(i) = a(i) + b(i) 920 continue return end c%4.11 subroutine s4111(a,b,ip,n) c c indirect addressing c indirect addressing on lhs c integer n,ip(*) real a(*),b(*) do 560 i = 1,n a(ip(i)) = b(i) 560 continue return end c%4.11 subroutine s4112(a,b,ip,n) c c indirect addressing c indirect addressing on rhs c integer n,ip(*) real a(*),b(*) do 550 i = 1,n a(i) = b(ip(i)) 550 continue return end c%4.11 subroutine s4113(a,b,ip,n) c c indirect addressing c indirect addressing on rhs and lhs c integer n,ip(*) real a(*),b(*) do 570 i = 1,n a(ip(i)) = b(ip(i)) 570 continue return end c%4.11 subroutine s4114(xx,yy,n) c c indirect addressing c mixed variable and constant lb, ub, and stride c@ n1,n2 are not defined. c@ l is used before being defined c dimension xx(100), yy(100), zz(100), l(100) xx(1) = 5. j = 5 do 50 i=n1,n2 k=l(i) xx(i)=yy(i)*zz(n2-k+1) k = k + j 50 continue write (unit=6,fmt=11) j, xx(n) 11 format (i6,e12.6) return end c%4.11 subroutine s4115(xx,yy,n) c c indirect addressing c all variable - lb, ub, and stride c@ will always vectorize since I is induction variable and c@ XX appears only on left hand side. c@ n1,n2,n3 are not defined. c@ l is used before being defined c dimension xx(100), yy(100), zz(100), l(100) xx(1) = 7. do 70 i=n1,n2,n3 k=l(i) xx(i)=yy(2*i+1)*zz(k+i/n3+n1) 70 continue write (unit=6,fmt=11) xx(n) 11 format (i6,e12.6) return end c%4.11 subroutine s4116(n,d) c c indirect addressing c This example test partial vectorization for the creation of c vector of indices. The best results are partial vectorization c with the net effect being the loops being compiled as: c c temp(1)=1 c do i= 2,n c temp(i) = temp(i-1)*2 (actually temp(i-1)+temp(i-1)) c enddo c do i= 1,n c d(temp(i)) =1 c enddo c dimension d(*) j =1 do 40, i= 1,n j = j*2 d(j) = 1 40 continue return end c%4.11 subroutine s4117(a,b,c,n,m) c c indirect addressing c seq function c dimension a(1000),b(1000),c(1000) do 1200 i = 1,m a(i) = b(i) + c(i/2) 1200 continue return end .