54 real(dp),
parameter::
zero=0
56 interface avco;
module procedure avco, davco, tavco; end
interface 57 interface dfco;
module procedure dfco, ddfco, tdfco; end
interface 58 interface dfco2;
module procedure dfco2, ddfco2, tdfco2; end
interface 64 interface ldub;
module procedure ldub, dldub; end
interface 65 interface ldltb;
module procedure ldltb, dldltb; end
interface 66 interface l1ubb;
module procedure l1ubb, dl1ubb; end
interface 67 interface l1ueb;
module procedure l1ueb, dl1ueb; end
interface 69 interface udlb;
module procedure udlb, dudlb; end
interface 83 interface wrtb;
module procedure wrtb; end
interface 99 subroutine avco(na,nb,za,zb,z0,a,b)
103 integer(spi),
intent(in ):: na,nb
104 real(sp),
intent(in ):: za(na),zb(nb),z0
105 real(sp),
intent(out):: a(na),b(nb)
107 integer(spi) :: na1,nab,i
108 real(sp),
dimension(na+nb,na+nb):: w
109 real(sp),
dimension(na) :: za0,pa
110 real(sp),
dimension(nb) :: zb0,pb
111 real(sp),
dimension(na+nb) :: ab
117 w(1,1:na)=
u1; ab(1)=
u1 118 do i=2,nab; w(i,1:na)=pa; pa=pa*za0; w(i,na1:nab)=pb; pb=pb*zb0;
enddo 120 a=ab(1:na); b=ab(na1:nab)
133 subroutine davco(na,nb,za,zb,z0,a,b)
137 integer(spi),
intent(IN ):: na,nb
138 real(dp),
intent(IN ):: za(na),zb(nb),z0
139 real(dp),
intent(OUT):: a(na),b(nb)
141 integer(spi) :: na1,nab,i
142 real(dp),
dimension(na+nb,na+nb):: w
143 real(dp),
dimension(na) :: za0,pa
144 real(dp),
dimension(nb) :: zb0,pb
145 real(dp),
dimension(na+nb) :: ab
151 w(1,1:na)=
u1; ab(1)=
u1 152 do i=2,nab; w(i,1:na)=pa; pa=pa*za0; w(i,na1:nab)=pb; pb=pb*zb0;
enddo 154 a=ab(1:na); b=ab(na1:nab)
164 subroutine tavco(xa,xb,a,b)
166 real(dp),
dimension(:),
intent(IN ):: xa,xb
167 real(dp),
dimension(:),
intent(OUT):: a,b
171 na=
size(xa);
if(na /=
size(a))stop
'In tavco; sizes of a and xa different' 172 nb=
size(xb);
if(nb /=
size(b))stop
'In tavco; sizes of b and xb different' 173 call davco(na,nb,xa,xb,
zero,a,b)
189 subroutine dfco(na,nb,za,zb,z0,a,b)
190 use pietc_s,
only: u0,u1
193 integer(spi),
intent(IN ) :: na,nb
194 real(sp),
intent(IN ) :: za(na),zb(nb),z0
195 real(sp),
intent(OUT) :: a(na),b(nb)
197 integer(spi):: na1,nab,i
198 real(sp),
dimension(na+nb,na+nb):: w
199 real(sp),
dimension(na) :: za0,pa
200 real(sp),
dimension(nb) :: zb0,pb
201 real(sp),
dimension(na+nb) :: ab
207 w(1,1:na)=u1; ab(1)=u1
208 do i=3,nab; w(i,1:na) =pa*(i-2); pa=pa*za0;
enddo 209 do i=2,nab; w(i,na1:nab)=pb; pb=pb*zb0;
enddo 211 a=ab(1:na); b=ab(na1:nab)
224 subroutine ddfco(na,nb,za,zb,z0,a,b)
228 integer(spi),
intent(in) :: na,nb
229 real(dp),
intent(in) :: za(na),zb(nb),z0
230 real(dp),
intent(out):: a(na),b(nb)
232 integer(spi) :: na1,nab,i
233 real(dp),
dimension(na+nb,na+nb):: w
234 real(dp),
dimension(na) :: za0,pa
235 real(dp),
dimension(nb) :: zb0,pb
236 real(dp),
dimension(na+nb) :: ab
242 w(1,1:na)=
u1; ab(1)=
u1 243 do i=3,nab; w(i,1:na) =pa*(i-2); pa=pa*za0;
enddo 244 do i=2,nab; w(i,na1:nab)=pb; pb=pb*zb0;
enddo 246 a=ab(1:na); b=ab(na1:nab)
257 subroutine tdfco(xa,xb,a,b)
259 real(dp),
dimension(:),
intent(IN ):: xa,xb
260 real(dp),
dimension(:),
intent(OUT):: a,b
264 na=
size(xa);
if(na /=
size(a))stop
'In tdfco; sizes of a and xa different' 265 nb=
size(xb);
if(nb /=
size(b))stop
'In tdfco; sizes of b and xb different' 266 call ddfco(na,nb,xa,xb,
zero,a,b)
282 subroutine dfco2(na,nb,za,zb,z0,a,b)
283 use pietc_s,
only: u0,u1
286 integer(spi),
intent(IN ):: na,nb
287 real(sp),
intent(IN ):: za(na),zb(nb),z0
288 real(sp),
intent(OUT):: a(na),b(nb)
290 integer(spi) :: na1,nab,i
291 real(sp),
dimension(na+nb,na+nb):: w
292 real(sp),
dimension(na) :: za0,pa
293 real(sp),
dimension(nb) :: zb0,pb
294 real(sp),
dimension(na+nb) :: ab
300 w(1,1:na)=u1; ab(1)=u1
301 do i=4,nab; w(i,1:na) =pa*(i-2)*(i-3); pa=pa*za0;
enddo 302 do i=2,nab; w(i,na1:nab)=pb; pb=pb*zb0;
enddo 304 a=ab(1:na); b=ab(na1:nab)
317 subroutine ddfco2(na,nb,za,zb,z0,a,b)
321 integer(spi),
intent(IN ) :: na,nb
322 real(dp),
intent(IN ) :: za(na),zb(nb),z0
323 real(dp),
intent(OUT) :: a(na),b(nb)
325 integer(spi) :: na1,nab,i
326 real(dp),
dimension(na+nb,na+nb):: w
327 real(dp),
dimension(na) :: za0,pa
328 real(dp),
dimension(nb) :: zb0,pb
329 real(dp),
dimension(na+nb) :: ab
335 w(1,1:na)=
u1; ab(1)=
u1 336 do i=4,nab; w(i,1:na) =pa*(i-2)*(i-3); pa=pa*za0;
enddo 337 do i=2,nab; w(i,na1:nab)=pb; pb=pb*zb0;
enddo 339 a=ab(1:na); b=ab(na1:nab)
340 end subroutine ddfco2
349 subroutine tdfco2(xa,xb,a,b)
351 real(dp),
dimension(:),
intent(IN ):: xa,xb
352 real(dp),
dimension(:),
intent(OUT):: a,b
356 na=
size(xa);
if(na /=
size(a))stop
'In tdfco2; sizes of a and xa different' 357 nb=
size(xb);
if(nb /=
size(b))stop
'In tdfco2; sizes of b and xb different' 358 call ddfco2(na,nb,xa,xb,
zero,a,b)
359 end subroutine tdfco2
369 pure subroutine clib(m1,m2,mah1,mah2,a)
370 use pietc_s,
only: u0
372 integer(spi),
intent(IN ) :: m1, m2, mah1, mah2
373 real(sp),
intent(INOUT) :: a(m1,-mah1:mah2)
375 do j=1,mah1; a(1:min(m1,j),-j)=u0;
enddo 376 do j=m2-m1+1,mah2; a(max(1,m2-j+1):m1,j)=u0;
enddo 387 pure subroutine clib_d(m1,m2,mah1,mah2,a)
390 integer(spi),
intent(IN ) :: m1, m2, mah1, mah2
391 real(dp),
intent(INOUT) :: a(m1,-mah1:mah2)
393 do j=1,mah1; a(1:min(m1,j),-j)=
u0;
enddo 394 do j=m2-m1+1,mah2; a(max(1,m2-j+1):m1,j)=
u0;
enddo 405 pure subroutine clib_c(m1,m2,mah1,mah2,a)
408 integer(spi),
intent(IN ) :: m1, m2, mah1, mah2
409 complex(dpc),
intent(INOUT) :: a(m1,-mah1:mah2)
411 do j=1,mah1; a(1:min(m1,j),-j)=
c0;
enddo 412 do j=m2-m1+1,mah2; a(max(1,m2-j+1):m1,j)=
c0;
enddo 425 subroutine cad1b(m1,mah1,mah2,mirror2,a)
426 use pietc_s,
only: u0
428 integer(spi),
intent(IN ):: m1,mah1,mah2,mirror2
429 real(sp),
intent(INOUT):: a(0:m1-1,-mah1:mah2)
431 integer(spi):: i,i2,jm,jp,jpmax
433 if(mirror2+mah1 > mah2)stop
'In CAD1B; mah2 insufficient' 434 do i=0,m1-1; i2=i*2; jpmax=mirror2+mah1-i2;
if(jpmax <= -mah1)
exit 435 do jm=-mah1,mah2; jp=mirror2-jm-i2;
if(jp <= jm)
exit 436 a(i,jp)=a(i,jp)+a(i,jm)
450 subroutine csb1b(m1,mah1,mah2,mirror2,a)
451 use pietc_s,
only: u0
453 integer(spi),
intent(IN ):: m1,mah1,mah2,mirror2
454 real(sp),
intent(INOUT):: a(0:m1-1,-mah1:mah2)
456 integer(spi):: i,i2,jm,jp,jpmax
458 if(mirror2+mah1 > mah2)stop
'In CSB1B; mah2 insufficient' 459 do i=0,m1-1; i2=i*2; jpmax=mirror2+mah1-i2;
if(jpmax < -mah1)
exit 460 do jm=-mah1,mah2; jp=mirror2-jm-i2;
if(jp < jm)
exit 461 a(i,jp)=a(i,jp)-a(i,jm)
476 subroutine cad2b(m1,m2,mah1,mah2,mirror2,a)
477 use pietc_s,
only: u0
479 integer(spi),
intent(IN ):: m1,m2,mah1,mah2,mirror2
480 real(sp),
intent(INOUT):: a(1-m1:0,m1-m2-mah1:m1-m2+mah2)
482 integer(spi):: i,i2,jm,jp,jmmin,nah1,nah2
484 nah1=mah1+m2-m1; nah2=mah2+m1-m2
485 if(mirror2-nah1 > -nah2)stop
'In CAD2B; mah1 insufficient' 486 do i=0,1-m1,-1; i2=i*2; jmmin=mirror2-nah2-i2;
if(jmmin >= nah2)
exit 487 do jp=nah2,nah1,-1; jm=mirror2-jp-i2;
if(jm >= jp)
exit 488 a(i,jm)=a(i,jm)+a(i,jp)
503 subroutine csb2b(m1,m2,mah1,mah2,mirror2,a)
504 use pietc_s,
only: u0
506 integer(spi),
intent(IN ):: m1,m2,mah1,mah2,mirror2
507 real(sp),
intent(INOUT):: a(1-m1:0,m1-m2-mah1:m1-m2+mah2)
509 integer(spi):: i,i2,jm,jp,jmmin,nah1,nah2
511 nah1=mah1+m2-m1; nah2=mah2+m1-m2
512 if(mirror2-nah1 > -nah2)stop
'In CSB2B; mah1 insufficient' 513 do i=0,1-m1,-1; i2=i*2; jmmin=mirror2-nah2-i2;
if(jmmin > nah2)
exit 514 do jp=nah2,nah1,-1; jm=mirror2-jp-i2;
if(jm > jp)
exit 515 a(i,jm)=a(i,jm)-a(i,jp)
531 subroutine ldub(m,mah1,mah2,a)
532 use pietc_s,
only: u0,u1
534 integer(spi),
intent(IN ):: m,mah1, mah2
535 real(sp),
intent(INOUT):: a(m,-mah1:mah2)
537 integer(spi):: j, imost, jmost, jp, i
538 real(sp) :: ajj, ajji, aij
546 print
'(" Failure in LDUB:"/" Matrix requires pivoting or is singular")' 554 a(i,jp-i:jmost-i)=a(i,jp-i:jmost-i)-aij*a(j,1:jmost-j)
556 a(j,1:jmost-j)=ajji*a(j,1:jmost-j)
570 subroutine dldub(m,mah1,mah2,a)
573 integer(spi),
intent(IN ):: m,mah1, mah2
574 real(dp),
intent(INOUT):: a(m,-mah1:mah2)
576 integer(spi):: j, imost, jmost, jp, i
577 real(dp) :: ajj, ajji, aij
585 print
'(" Fails in LDUB_d:"/" Matrix requires pivoting or is singular")' 593 a(i,jp-i:jmost-i)=a(i,jp-i:jmost-i)-aij*a(j,1:jmost-j)
595 a(j,1:jmost-j)=ajji*a(j,1:jmost-j)
608 subroutine ldltb(m,mah1,a)
609 use pietc_s,
only: u0,u1
610 integer(spi),
intent(IN ):: m,mah1
611 real(sp),
intent(INOUT):: a(m,-mah1:0)
613 integer(spi):: j, imost, jp, i,k
614 real(sp) :: ajj, ajji, aij
621 print
'(" Fails in LDLTB:"/" Matrix requires pivoting or is singular")' 630 a(i,k-i)=a(i,k-i)-aij*a(k,j-k)
644 subroutine dldltb(m,mah1,a)
646 integer(spi),
intent(IN ) :: m,mah1
647 real(dp),
intent(INOUT) :: a(m,-mah1:0)
649 integer(spi):: j, imost, jp, i,k
650 real(dp) :: ajj, ajji, aij
657 print
'(" Fails in LDLTB_d:"/" Matrix requires pivoting or is singular")' 666 a(i,k-i)=a(i,k-i)-aij*a(k,j-k)
670 end subroutine dldltb
682 subroutine udlb(m,mah1,mah2,a)
684 integer(spi),
intent(IN ) :: m,mah1,mah2
685 real(sp),
dimension(m,-mah1:mah2),
intent(INOUT) :: a(m,-mah1:mah2)
687 real(sp),
dimension(m,-mah2:mah1):: at
689 at=a(m:1:-1,mah2:-mah1:-1);
call ldub(m,mah2,mah1,at)
690 a=at(m:1:-1,mah1:-mah2:-1)
703 subroutine dudlb(m,mah1,mah2,a)
705 integer(spi),
intent(IN ) :: m,mah1,mah2
706 real(dp),
dimension(m,-mah1:mah2),
intent(INOUT) :: a(m,-mah1:mah2)
708 real(dp),
dimension(m,-mah2:mah1):: at
710 at=a(m:1:-1,mah2:-mah1:-1);
call dldub(m,mah2,mah1,at)
711 a=at(m:1:-1,mah1:-mah2:-1)
732 subroutine l1ubb(m,mah1,mah2,mbh1,mbh2,a,b)
733 use pietc_s,
only: u0,u1
735 integer(spi),
intent(IN ) :: m,mah1, mah2, mbh1, mbh2
736 real(sp),
intent(INOUT) :: a(m,-mah1:mah2), b(m,-mbh1:mbh2)
738 integer(spi):: j, imost, jmost, jleast, jp, i
739 real(sp) :: ajj, ajji, aij
747 if(ajj == u0)stop
'In L1UBB; zero element found in diagonal factor' 749 a(j,jleast-j:jmost-j) = ajji * a(j,jleast-j:jmost-j)
752 a(i,jp-i:jmost-i) = a(i,jp-i:jmost-i) - aij*a(j,jp-j:jmost-j)
755 b(j,-mbh1:mbh2) = ajji * b(j,-mbh1:mbh2)
770 subroutine dl1ubb(m,mah1,mah2,mbh1,mbh2,a,b)
773 integer(spi),
intent(IN ) :: m,mah1, mah2, mbh1, mbh2
774 real(dp),
intent(INOUT) :: a(m,-mah1:mah2), b(m,-mbh1:mbh2)
776 integer(spi):: j, imost, jmost, jleast, jp, i
777 real(dp) :: ajj, ajji, aij
785 if(ajj ==
u0)stop
'In L1UBB_d; zero element found in diagonal factor' 787 a(j,jleast-j:jmost-j) = ajji * a(j,jleast-j:jmost-j)
790 a(i,jp-i:jmost-i) = a(i,jp-i:jmost-i) - aij*a(j,jp-j:jmost-j)
793 b(j,-mbh1:mbh2) = ajji * b(j,-mbh1:mbh2)
795 end subroutine dl1ubb
817 subroutine l1ueb(m,mah1,mah2,mbh1,mbh2,a,b)
818 use pietc_s,
only: u0,u1
820 integer(spi),
intent(IN ) :: m,mah1, mah2, mbh1, mbh2
821 real(sp),
intent(INOUT) :: a(0:m,-mah1:mah2), b(m,-mbh1:mbh2)
823 integer(spi):: j, imost, jmost, jleast, jp, i
824 real(sp) :: ajj, ajji, aij
832 if(ajj == u0)stop
'In L1UEB; zero element found in diagonal factor' 834 a(j,jleast-j:jmost-j) = ajji * a(j,jleast-j:jmost-j)
837 a(i,jp-i:jmost-i) = a(i,jp-i:jmost-i) - aij*a(j,jp-j:jmost-j)
840 b(j,-mbh1:mbh2) = ajji * b(j,-mbh1:mbh2)
854 subroutine dl1ueb(m,mah1,mah2,mbh1,mbh2,a,b)
857 integer(spi),
intent(IN ):: m,mah1, mah2, mbh1, mbh2
858 real(dp),
intent(INOUT):: a(0:,-mah1:), b(:,-mbh1:)
860 integer(spi):: j, imost, jmost, jleast, jp, i
861 real(dp) :: ajj, ajji, aij
869 if(ajj ==
u0)stop
'In L1UEB_D; zero element found in diagonal factor' 871 a(j,jleast-j:jmost-j) = ajji * a(j,jleast-j:jmost-j)
874 a(i,jp-i:jmost-i) = a(i,jp-i:jmost-i) - aij*a(j,jp-j:jmost-j)
877 b(j,-mbh1:mbh2) = ajji * b(j,-mbh1:mbh2)
879 end subroutine dl1ueb
891 subroutine udlbv(m,mah1,mah2,a,v)
893 integer(spi),
intent(IN ):: m, mah1, mah2
894 real(sp),
intent(IN ):: a(m,-mah1:mah2)
895 real(sp),
intent(INOUT):: v(m)
902 do i=j+1,min(m,j+mah1); v(i)=v(i)-a(i,j-i)*vj; enddo; v(j)=a(j,0)*vj
906 do i=max(1,j-mah2),j-1; v(i)=v(i)-a(i,j-i)*vj;
enddo 920 subroutine dudlbv(m,mah1,mah2,a,v)
922 integer(spi),
intent(IN ) :: m, mah1, mah2
923 real(dp),
intent(IN ) :: a(m,-mah1:mah2)
924 real(dp),
intent(INOUT) :: v(m)
931 do i=j+1,min(m,j+mah1); v(i)=v(i)-a(i,j-i)*vj; enddo; v(j)=a(j,0)*vj
935 do i=max(1,j-mah2),j-1; v(i)=v(i)-a(i,j-i)*vj;
enddo 949 subroutine ltdlbv(m,mah1,a,v)
951 integer(spi),
intent(IN ) :: m, mah1
952 real(sp),
intent(IN ) :: a(m,-mah1:0)
953 real(sp),
intent(INOUT) :: v(m)
960 do i=j+1,min(m,j+mah1); v(i)=v(i)-a(i,j-i)*vj; enddo; v(j)=a(j,0)*vj
964 do i=max(1,j-mah1),j-1; v(i)=v(i)-a(j,i-j)*vj;
enddo 981 integer(spi),
intent(IN ) :: m, mah1
982 real(dp),
intent(IN ) :: a(m,-mah1:0)
983 real(dp),
intent(INOUT) :: v(m)
990 do i=j+1,min(m,j+mah1); v(i)=v(i)-a(i,j-i)*vj; enddo; v(j)=a(j,0)*vj
994 do i=max(1,j-mah1),j-1; v(i)=v(i)-a(j,i-j)*vj;
enddo 1010 subroutine udlbx(mx,mah1,mah2,my,a,v)
1012 integer(spi),
intent(IN ) :: mx, mah1, mah2, my
1013 real(sp),
intent(IN ) :: a(mx,-mah1:mah2)
1014 real(sp),
intent(INOUT) :: v(mx,my)
1016 integer(spi):: jx, ix
1019 do ix=jx+1,min(mx,jx+mah1); v(ix,:) = v(ix,:) - a(ix,jx-ix)*v(jx,:);
enddo 1020 v(jx,:) = a(jx,0) * v(jx,:)
1023 do ix=max(1,jx-mah2),jx-1; v(ix,:) = v(ix,:) - a(ix,jx-ix)*v(jx,:);
enddo 1025 end subroutine udlbx 1039 subroutine udlby(my,mah1,mah2,mx,a,v)
1041 integer(spi),
intent(IN ) :: my, mah1, mah2, mx
1042 real(sp),
intent(IN ) :: a(my,-mah1:mah2)
1043 real(sp),
intent(INOUT) :: v(mx,my)
1045 integer(spi):: iy, jy
1048 do iy=jy+1,min(my,jy+mah1); v(:,iy) = v(:,iy)-a(iy,jy-iy)*v(:,jy);
enddo 1049 v(:,jy)=a(jy,0)*v(:,jy)
1052 do iy=max(1,jy-mah2),jy-1; v(:,iy)=v(:,iy)-a(iy,jy-iy)*v(:,jy);
enddo 1054 end subroutine udlby 1066 subroutine udlvb(m,mah1,mah2,v,a)
1068 integer(spi),
intent(IN ) :: m, mah1, mah2
1069 real(sp),
intent(IN ) :: a(m,-mah1:mah2)
1070 real(sp),
intent(INOUT) :: v(m)
1077 do j=i+1,min(m,i+mah2); v(j)=v(j)-vi*a(i,j-i);
enddo 1082 do j=max(1,i-mah1),i-1; v(j)=v(j)-vi*a(i,j-i);
enddo 1084 end subroutine udlvb 1098 subroutine udlxb(mx,mah1,mah2,my,v,a)
1100 integer(spi),
intent(IN ) :: mx, mah1, mah2, my
1101 real(sp),
intent(IN ) :: a(mx,-mah1:mah2)
1102 real(sp),
intent(INOUT) :: v(mx,my)
1104 integer(spi):: ix, jx
1107 do jx=ix+1,min(mx,ix+mah2); v(jx,:)=v(jx,:)-v(ix,:)*a(ix,jx-ix);
enddo 1108 v(ix,:)=v(ix,:)*a(ix,0)
1111 do jx=max(1,ix-mah1),ix-1; v(jx,:)=v(jx,:)-v(ix,:)*a(ix,jx-ix);
enddo 1113 end subroutine udlxb 1127 subroutine udlyb(my,mah1,mah2,mx,v,a)
1129 integer(spi),
intent(IN ) :: my, mah1, mah2, mx
1130 real(sp),
intent(IN ) :: a(my,-mah1:mah2)
1131 real(sp),
intent(INOUT) :: v(mx,my)
1133 integer(spi):: iy, jy
1136 do jy=iy+1,min(my,iy+mah2); v(:,jy)=v(:,jy)-v(:,iy)*a(iy,jy-iy);
enddo 1137 v(:,iy)=v(:,iy)*a(iy,0)
1140 do jy=max(1,iy-mah1),iy-1; v(:,jy)=v(:,jy)-v(:,iy)*a(iy,jy-iy);
enddo 1142 end subroutine udlyb 1154 subroutine u1lbv(m,mah1,mah2,a,v)
1156 integer(spi),
intent(IN ) :: m, mah1, mah2
1157 real(sp),
intent(IN ) :: a(m,-mah1:mah2)
1158 real(sp),
intent(INOUT) :: v(m)
1165 do i=j+1,min(m,j+mah1); v(i)=v(i)-a(i,j-i)*vj;
enddo 1169 do i=max(1,j-mah2),j-1; v(i)=v(i)-a(i,j-i)*vj;
enddo 1171 end subroutine u1lbv 1185 subroutine u1lbx(mx,mah1,mah2,my,a,v)
1187 integer(spi),
intent(IN ) :: mx, mah1, mah2, my
1188 real(sp),
intent(IN ) :: a(mx,-mah1:mah2)
1189 real(sp),
intent(INOUT) :: v(mx,my)
1191 integer(spi):: ix, jx
1194 do ix=jx+1,min(mx,jx+mah1); v(ix,:)=v(ix,:)-a(ix,jx-ix)*v(jx,:);
enddo 1197 do ix=max(1,jx-mah2),jx-1; v(ix,:)=v(ix,:)-a(ix,jx-ix)*v(jx,:);
enddo 1199 end subroutine u1lbx 1213 subroutine u1lby(my,mah1,mah2,mx,a,v)
1215 integer(spi),
intent(IN ) :: my, mah1, mah2, mx
1216 real(sp),
intent(IN ) :: a(my,-mah1:mah2)
1217 real(sp),
intent(INOUT) :: v(mx,my)
1219 integer(spi):: iy, jy
1222 do iy=jy+1,min(my,jy+mah1); v(:,iy)=v(:,iy)-a(iy,jy-iy)*v(:,jy);
enddo 1225 do iy=max(1,jy-mah2),jy-1; v(:,iy)=v(:,iy)-a(iy,jy-iy)*v(:,jy);
enddo 1227 end subroutine u1lby 1239 subroutine u1lvb(m,mah1,mah2,v,a)
1241 integer(spi),
intent(IN ) :: m, mah1, mah2
1242 real(sp),
intent(IN ) :: a(m,-mah1:mah2)
1243 real(sp),
intent(INOUT) :: v(m)
1250 do j=i+1,min(m,i+mah2); v(j)=v(j)-vi*a(i,j-i);
enddo 1254 do j=max(1,i-mah1),i-1; v(j)=v(j)-vi*a(i,j-i);
enddo 1256 end subroutine u1lvb 1270 subroutine u1lxb(mx,mah1,mah2,my,v,a)
1272 integer(spi),
intent(IN ) :: mx, mah1, mah2, my
1273 real(sp),
intent(IN ) :: a(mx,-mah1:mah2)
1274 real(sp),
intent(INOUT) :: v(mx,my)
1276 integer(spi):: ix, jx
1279 do jx=ix+1,min(mx,ix+mah2); v(jx,:)=v(jx,:)-v(ix,:)*a(ix,jx-ix);
enddo 1282 do jx=max(1,ix-mah1),ix-1; v(jx,:)=v(jx,:)-v(ix,:)*a(ix,jx-ix);
enddo 1284 end subroutine u1lxb 1298 subroutine u1lyb(my,mah1,mah2,mx,v,a)
1300 integer(spi),
intent(IN ) :: my, mah1, mah2, mx
1301 real(sp),
intent(IN ) :: a(my,-mah1:mah2)
1302 real(sp),
intent(INOUT) :: v(mx,my)
1304 integer(spi):: iy, jy
1307 do jy=iy+1,min(my,iy+mah2); v(:,jy)=v(:,jy)-v(:,iy)*a(iy,jy-iy);
enddo 1310 do jy=max(1,iy-mah1),iy-1; v(:,jy)=v(:,jy)-v(:,iy)*a(iy,jy-iy);
enddo 1312 end subroutine u1lyb 1322 subroutine linbv(m,mah1,mah2,a,v)
1324 integer(spi),
intent(IN ) :: m, mah1, mah2
1325 real(sp),
intent(INOUT) :: a(m,-mah1:mah2), v(m)
1327 call ldub(m,mah1,mah2,a)
1328 call udlbv(m,mah1,mah2,a,v)
1329 end subroutine linbv 1340 subroutine wrtb(m1,m2,mah1,mah2,a)
1342 integer(spi),
intent(IN) :: m1, m2, mah1, mah2
1343 real(sp),
intent(IN) :: a(m1,-mah1:mah2)
1345 integer(spi):: i1, i2, i, j1, j2, j, nj1
1349 print
'(7x,6(i2,10x))',(j,j=-mah1,mah2)
1354 if(nj1==0)print
'(1x,i3,6(1x,e12.5))', i,(a(i,j),j=j1,j2)
1355 if(nj1==1)print
'(1x,i3,12x,5(1x,e12.5))',i,(a(i,j),j=j1,j2)
1356 if(nj1==2)print
'(1x,i3,24x,4(1x,e12.5))',i,(a(i,j),j=j1,j2)
1357 if(nj1==3)print
'(1x,i3,36x,3(1x,e12.5))',i,(a(i,j),j=j1,j2)
1358 if(nj1==4)print
'(1x,i3,48x,2(1x,e12.5))',i,(a(i,j),j=j1,j2)
1359 if(nj1==5)print
'(1x,i3,60x,1(1x,e12.5))',i,(a(i,j),j=j1,j2)
integer, parameter sp
Single precision real kind.
integer, parameter dp
Double precision real kind.
subroutine dudlbv(m, mah1, mah2, a, v)
Back-substitution step of linear inversion involving banded LDU factored matrix and input vector...
Standard integer, real, and complex single and double precision kinds.
real(dp), parameter u0
zero
pure subroutine clib_d(m1, m2, mah1, mah2, a)
Clip (set to zero) the unused values in a banded matrix representation.
real(dp), parameter zero
Double precision real zero.
subroutine dltdlbv(m, mah1, a, v)
Like udlbv, except assuming a is the ltdl decomposition of a SYMMETRIC banded matrix, with only the non-upper part provided (to avoid redundancy)
complex(dpc), parameter c0
complex zero
real(dp), parameter u1
one
pure subroutine clib_c(m1, m2, mah1, mah2, a)
Clip (set to zero) the unused values in a banded matrix representation.
Some of the commonly used constants (pi etc) mainly for double-precision subroutines.
integer, parameter spi
Single precision integer kind.
integer, parameter dpc
Double precision real kind.
Routines dealing with the operations of banded matrices.