105interface ctoz;
module procedure ctoz;
end interface
120real(
sp),
dimension(:),
intent(in):: a
122s=sqrt(dot_product(a,a))
132real(
dp),
dimension(:),
intent(in):: a
134s=sqrt(dot_product(a,a))
145real(
sp),
dimension(:),
intent(IN):: a
146real(
sp),
dimension(size(a)) :: b
148s=
absv_s(a);
if(s==u0)then; b=u0;else;b=a/s;
endif
159real(
dp),
dimension(:),
intent(IN):: a
160real(
dp),
dimension(size(a)) :: b
162s=
absv_d(a);
if(s==
u0)then; b=
u0;else;b=a/s;
endif
173real(
sp),
dimension(:),
intent(in):: u,a
174real(
sp),
dimension(size(u)) :: b
177s=dot_product(u,a); b=a-u*s
188real(
dp),
dimension(:),
intent(in):: u,a
189real(
dp),
dimension(size(u)) :: b
192s=dot_product(u,a); b=a-u*s
203real(
sp),
dimension(3),
intent(in):: a,b
204real(
sp),
dimension(3) :: c
205c(1)=a(2)*b(3)-a(3)*b(2); c(2)=a(3)*b(1)-a(1)*b(3); c(3)=a(1)*b(2)-a(2)*b(1)
216real(
dp),
dimension(3),
intent(in):: a,b
217real(
dp),
dimension(3) :: c
218c(1)=a(2)*b(3)-a(3)*b(2); c(2)=a(3)*b(1)-a(1)*b(3); c(3)=a(1)*b(2)-a(2)*b(1)
233real(
sp),
dimension(4),
intent(in ):: u,v,w
234real(
sp),
dimension(4) :: x
235real(
sp):: uv12,uv13,uv14,uv23,uv24,uv34
236uv12=u(1)*v(2)-u(2)*v(1); uv13=u(1)*v(3)-u(3)*v(1); uv14=u(1)*v(4)-u(4)*v(1)
237 uv23=u(2)*v(3)-u(3)*v(2); uv24=u(2)*v(4)-u(4)*v(2)
238 uv34=u(3)*v(4)-u(4)*v(3)
239x(1)=-uv23*w(4)+uv24*w(3)-uv34*w(2)
240x(2)= uv13*w(4)-uv14*w(3) +uv34*w(1)
241x(3)=-uv12*w(4) +uv14*w(2)-uv24*w(1)
242x(4)= uv12*w(3)-uv13*w(2)+uv23*w(1)
254real(
dp),
dimension(4),
intent(in ):: u,v,w
255real(
dp),
dimension(4) :: x
256real(
dp):: uv12,uv13,uv14,uv23,uv24,uv34
257uv12=u(1)*v(2)-u(2)*v(1); uv13=u(1)*v(3)-u(3)*v(1); uv14=u(1)*v(4)-u(4)*v(1)
258 uv23=u(2)*v(3)-u(3)*v(2); uv24=u(2)*v(4)-u(4)*v(2)
259 uv34=u(3)*v(4)-u(4)*v(3)
260x(1)=-uv23*w(4)+uv24*w(3)-uv34*w(2)
261x(2)= uv13*w(4)-uv14*w(3) +uv34*w(1)
262x(3)=-uv12*w(4) +uv14*w(2)-uv24*w(1)
263x(4)= uv12*w(3)-uv13*w(2)+uv23*w(1)
274real(
sp),
dimension(:),
intent(in ):: a
275real(
sp),
dimension(:),
intent(in ):: b
276real(
sp),
DIMENSION(size(a),size(b)):: c
279do i=1,nb; c(:,i)=a*b(i);
enddo
290real(
dp),
dimension(:),
intent(in ):: a
291real(
dp),
dimension(:),
intent(in ):: b
292real(
dp),
dimension(size(a),size(b)):: c
295do i=1,nb; c(:,i)=a*b(i);
enddo
306integer(spi),
dimension(:),
intent(in ):: a
307integer(spi),
dimension(:),
intent(in ):: b
308integer(spi),
dimension(size(a),size(b)):: c
311do i=1,nb; c(:,i)=a*b(i);
enddo
323real(
sp),
dimension(3),
intent(IN ):: a,b,c
324real(
sp) :: tripleproduct
337real(
dp),
dimension(3),
intent(IN ):: a,b,c
338real(
dp) :: tripleproduct
350real(
sp),
dimension(:,:),
intent(IN ) :: a
352real(
sp),
dimension(size(a,1),size(a,1)):: b
353integer(spi) :: n,nrank
371real(
dp),
dimension(:,:),
intent(IN ) :: a
373real(
dp),
dimension(size(a,1),size(a,1)):: b
374integer(spi) :: n,nrank
391integer(spi),
dimension(:,:),
intent(IN ):: a
393real(
dp),
dimension(size(a,1),size(a,2)):: b
395b=a; bdet=
det(b); idet=nint(bdet)
406integer(dpi),
dimension(:,:),
intent(IN ):: a
408real(
dp),
dimension(size(a,1),size(a,2)) :: b
410b=a; bdet=
det(b); idet=nint(bdet)
422real(sp),
dimension(3),
intent(IN ):: a
423real(sp),
dimension(3,3) :: b
424b=u0;b(3,2)=a(1);b(1,3)=a(2);b(2,1)=a(3);b(2,3)=-a(1);b(3,1)=-a(2);b(1,2)=-a(3)
436real(dp),
dimension(3),
intent(IN ):: a
437real(dp),
dimension(3,3) :: b
438b=
u0;b(3,2)=a(1);b(1,3)=a(2);b(2,1)=a(3);b(2,3)=-a(1);b(3,1)=-a(2);b(1,2)=-a(3)
450real(sp),
dimension(3,3),
intent(IN ):: b
451real(sp),
dimension(3) :: a
452a(1)=(b(3,2)-b(2,3))*o2; a(2)=(b(1,3)-b(3,1))*o2; a(3)=(b(2,1)-b(1,2))*o2
464real(dp),
dimension(3,3),
intent(IN ):: b
465real(dp),
dimension(3) :: a
466a(1)=(b(3,2)-b(2,3))*
o2; a(2)=(b(1,3)-b(3,1))*
o2; a(3)=(b(2,1)-b(1,2))*
o2
478real(sp),
dimension(:),
intent(IN ) :: a
479real(sp),
dimension(size(a),size(a)):: b
482b=
u0;
do i=1,n; b(i,i)=a(i);
enddo
494real(dp),
dimension(:),
intent(IN ) :: a
495real(dp),
dimension(size(a),size(a)):: b
498b=
u0;
do i=1,n; b(i,i)=a(i);
enddo
509integer(spi),
dimension(:),
intent(IN ) :: a
510integer(spi),
dimension(size(a),size(a)):: b
513b=0;
do i=1,n; b(i,i)=a(i);
enddo
524real(sp),
dimension(:,:),
intent(IN ):: b
525real(sp),
dimension(size(b,1)) :: a
528do i=1,n; a(i)=b(i,i);
enddo
539real(dp),
dimension(:,:),
intent(IN ):: b
540real(dp),
dimension(size(b,1)) :: a
543do i=1,n; a(i)=b(i,i);
enddo
554integer(spi),
dimension(:,:),
intent(IN ):: b
555integer(spi),
dimension(size(b,1)) :: a
558do i=1,n; a(i)=b(i,i);
enddo
568real(sp),
dimension(:,:),
intent(IN ):: b
580real(dp),
dimension(:,:),
intent(IN ):: b
592integer(spi),
dimension(:,:),
intent(IN ):: b
604integer(spi),
intent(IN ) :: n
605integer(spi),
dimension(n,n):: a
607a=0;
do i=1,n; a(i,i)=1;
enddo
616integer(spi),
dimension(3,3):: a
618a=0;
do i=1,3; a(i,i)=1;
enddo
631real(sp),
intent(IN ):: sa,sb
636area=asin(sa*sb/(1+ca*cb))
649real(dp),
intent(IN ):: sa,sb
654area=asin(sa*sb/(1+ca*cb))
668use pietc_s,
only: zero=>u0
670real(sp),
dimension(3),
intent(IN ):: v1,v2,v3
672real(sp) :: s123,a1,a2,b,d1,d2,d3
673real(sp),
dimension(3) :: u0,u1,u2,u3,x,y
679d1=dot_product(u3-u2,u3-u2)
680d2=dot_product(u1-u3,u1-u3)
681d3=dot_product(u2-u1,u2-u1)
684if(d3<d1 .or. d3<d2)
call cyclic(u1,u2,u3,d1,d2,d3)
685if(d3<d1 .or. d3<d2)
call cyclic(u1,u2,u3,d1,d2,d3)
690a1=-dot_product(x,u1-u0); a2= dot_product(x,u2-u0)
706 real(sp),
dimension(3),
intent(INOUT):: u1,u2,u3
707 real(sp),
intent(INOUT):: d1,d2,d3
708 real(sp),
dimension(3) :: ut
710 dt=d1; d1=d2; d2=d3; d3=dt
711 ut=u1; u1=u2; u2=u3; u3=ut
725real(dp),
dimension(3),
intent(IN ):: v1,v2,v3
727real(dp) :: s123,a1,a2,b,d1,d2,d3
728real(dp),
dimension(3) ::
u0,
u1,
u2,
u3,x,y
745a1=-dot_product(x,
u1-
u0); a2= dot_product(x,
u2-
u0)
761 real(dp),
dimension(3),
intent(INOUT):: u1,u2,u3
762 real(dp),
intent(INOUT):: d1,d2,d3
763 real(dp),
dimension(3) :: ut
765 dt=d1; d1=d2; d2=d3; d3=dt
766 ut=u1; u1=u2; u2=u3; u3=ut
785function dtarea_s(rlat,drlata,drlona,drlatb,drlonb)
result(area)
786use pietc_s,
only: u0,u1
788real(sp),
intent(in ):: rlat,drlata,drlona,drlatb,drlonb
790real(sp),
dimension(2):: x2a,x2b,xb,yb
791real(sp) :: sb,ssb,cb,xa,sa,ca,sc,cc
792call dlltoxy(rlat,drlata,drlona,x2a)
793call dlltoxy(rlat,drlatb,drlonb,x2b)
794ssb=dot_product(x2b,x2b); sb=sqrt(ssb)
795if(sb==u0)then; area=u0; return;
endif
800xa=dot_product(xb,x2a)
801sa=dot_product(yb,x2a)
824function dtarea_d(rlat,drlata,drlona,drlatb,drlonb)
result(area)
827real(dp),
intent(in ):: rlat,drlata,drlona,drlatb,drlonb
829real(dp),
dimension(2):: x2a,x2b,xb,yb
830real(dp) :: sb,ssb,cb,xa,sa,ca,sc,cc
831call dlltoxy(rlat,drlata,drlona,x2a)
832call dlltoxy(rlat,drlatb,drlonb,x2b)
833ssb=dot_product(x2b,x2b); sb=sqrt(ssb)
834if(sb==
u0)then; area=
u0; return;
endif
839xa=dot_product(xb,x2a)
840sa=dot_product(yb,x2a)
868 (rlat,drlata,drlona,drlatb,drlonb,drlatc,drlonc) result(area)
870real(sp),
intent(in ):: rlat,drlata,drlona,drlatb,drlonb,drlatc,drlonc
872area=
sarea(rlat,drlata,drlona,drlatb,drlonb)&
873 -
sarea(rlat,drlatc,drlonc,drlatb,drlonb)
897 (rlat,drlata,drlona,drlatb,drlonb,drlatc,drlonc) result(area)
899real(dp),
intent(in ):: rlat,drlata,drlona,drlatb,drlonb,drlatc,drlonc
901area=
sarea(rlat,drlata,drlona,drlatb,drlonb)&
902 -
sarea(rlat,drlatc,drlonc,drlatb,drlonb)
918real(sp),
intent(in ):: rlat,drlat,drlon
919real(sp),
dimension(2),
intent(out):: x2
922x2=(/clata*sin(drlon),sin(drlat)+u2*sin(rlat)*clata*
hav(drlon)/)
938real(dp),
intent(in ):: rlat,drlat,drlon
939real(dp),
dimension(2),
intent(out):: x2
942x2=(/clata*sin(drlon),sin(drlat)+
u2*sin(rlat)*clata*
hav(drlon)/)
953real(sp),
intent(in ):: t
966real(dp),
intent(in )::
t
976use pietc_s,
only: u0,u1
978real(sp),
dimension(:),
intent(inout):: v
980s=
absv(v);
if(s==0)then; v=u0; v(1)=u1; else; v=v/s;
endif
990real(dp),
dimension(:),
intent(inout):: v
992s=
absv(v);
if(s==
u0)then; v=0; v(1)=
u1; else; v=v/s;
endif
1007use pietc_s,
only: u0,u1
1009real(sp),
dimension(:,:),
intent(IN ) :: as
1010real(sp),
dimension(:,:),
intent(OUT) :: b
1011integer(spi),
intent(OUT) :: nrank
1012real(sp),
intent(OUT) :: det
1013real(sp),
parameter :: crit=1.e-5_sp
1014real(sp),
dimension(size(as,1),size(as,2)):: a
1015real(sp),
dimension(size(as,2),size(as,1)):: ab
1016real(sp),
dimension(size(as,1)) :: tv,w
1017real(sp) :: val,s,vcrit
1018integer(spi) :: i,j,k,l,m,n
1019integer(spi),
dimension(2) :: ii
1022if(n/=
size(b,1) .or. n/=
size(b,2))stop
'In gram; incompatible dimensions'
1035 ab(k:m,k:n)=matmul( transpose(a(:,k:m)),b(:,k:n) )
1036 ii =maxloc( abs( ab(k:m,k:n)) )+k-1
1037 val=maxval( abs( ab(k:m,k:n)) )
1050 w(k:n)=matmul( transpose(b(:,k:n)),tv )
1051 b(:,k)=matmul(b(:,k:n),w(k:n) )
1052 s=dot_product(b(:,k),b(:,k))
1059 s=dot_product(b(:,l),b(:,j))
1080real(dp),
dimension(:,:),
intent(IN ) :: as
1081real(dp),
dimension(:,:),
intent(OUT) :: b
1082integer(spi),
intent(OUT) :: nrank
1083real(dp),
intent(OUT) :: det
1084real(dp),
parameter :: crit=1.e-9_dp
1085real(dp),
dimension(size(as,1),size(as,2)):: a
1086real(dp),
dimension(size(as,2),size(as,1)):: ab
1087real(dp),
dimension(size(as,1)) :: tv,w
1088real(dp) :: val,s,vcrit
1089integer(spi) :: i,j,k,l,m,n
1090integer(spi),
dimension(2) :: ii
1093if(n/=
size(b,1) .or. n/=
size(b,2))stop
'In gram; incompatible dimensions'
1106 ab(k:m,k:n)=matmul( transpose(a(:,k:m)),b(:,k:n) )
1107 ii =maxloc( abs( ab(k:m,k:n)) )+k-1
1108 val=maxval( abs( ab(k:m,k:n)) )
1121 w(k:n)=matmul( transpose(b(:,k:n)),tv )
1122 b(:,k)=matmul(b(:,k:n),w(k:n) )
1123 s=dot_product(b(:,k),b(:,k))
1130 s=dot_product(b(:,l),b(:,j))
1152real(dp),
dimension(:,:),
intent(IN ) :: as
1153real(dp),
dimension(:,:),
intent(OUT) :: b
1154integer(spi),
intent(OUT) :: nrank
1155integer(spi),
intent(out) :: detsign
1156real(dp),
intent(OUT) :: ldet
1157real(dp),
parameter :: crit=1.e-9_dp
1158real(dp),
dimension(size(as,1),size(as,2)):: a
1159real(dp),
dimension(size(as,2),size(as,1)):: ab
1160real(dp),
dimension(size(as,1)) :: tv,w
1161real(dp) :: val,s,vcrit
1162integer(spi) :: i,j,k,l,m,n
1163integer(spi),
dimension(2) :: ii
1167if(n/=
size(b,1) .or. n/=
size(b,2))stop
'In gram; incompatible dimensions'
1180 ab(k:m,k:n)=matmul( transpose(a(:,k:m)),b(:,k:n) )
1181 ii =maxloc( abs( ab(k:m,k:n)) )+k-1
1182 val=maxval( abs( ab(k:m,k:n)) )
1195 w(k:n)=matmul( transpose(b(:,k:n)),tv )
1196 b(:,k)=matmul(b(:,k:n),w(k:n) )
1197 s=dot_product(b(:,k),b(:,k))
1212 s=dot_product(b(:,l),b(:,j))
1227use pietc_s,
only: u0
1229real(sp),
dimension(:,:),
intent(INOUT) :: b
1230integer(spi),
intent( OUT) :: nrank
1231real(sp),
parameter :: crit=1.e-5_sp
1232real(sp) :: val,vcrit
1233integer(spi) :: j,k,n
1234n=
size(b,1);
if(n/=
size(b,2))stop
'In gram; matrix needs to be square'
1243 val=sqrt(dot_product(b(:,k),b(:,k)))
1251 b(:,j)=b(:,j)-b(:,k)*dot_product(b(:,k),b(:,j))
1265real(dp),
dimension(:,:),
intent(INOUT):: b
1266integer(spi),
intent( OUT):: nrank
1267real(dp),
parameter:: crit=1.e-9_dp
1268real(dp) :: val,vcrit
1269integer(spi) :: j,k,n
1270n=
size(b,1);
if(n/=
size(b,2))stop
'In gram; matrix needs to be square'
1279 val=sqrt(dot_product(b(:,k),b(:,k)))
1287 b(:,j)=b(:,j)-b(:,k)*dot_product(b(:,k),b(:,j))
1317integer(spi),
intent(IN ):: m,n
1318real(dp),
dimension(m,n),
intent(in ):: a
1319integer(spi),
dimension(n),
intent(out):: ipiv
1320real(dp),
dimension(m,n),
intent(out):: tt
1321real(dp),
dimension(n,n),
intent(out):: b
1322integer(spi),
intent(out):: rank
1323real(dp),
parameter :: eps=1.e-13_dp,epss=eps**2
1324real(dp),
dimension(m,n) :: aa
1325real(dp),
dimension(n) :: rowv
1326real(dp),
dimension(m) :: p
1327real(dp) :: maxp,nepss
1328integer(spi),
dimension(1):: jloc
1329integer(spi) :: i,ii,iii,j,maxi
1330if(m<n)stop
'In rowgram; this routines needs m>=n please'
1342 p(i)=dot_product(aa(i,:),aa(i,:))
1350 b(1:ii-1,:)=aa(1:ii-1,:)
1355 rowv(j)=maxval(abs(b(1:iii-1,j)))
1361 maxp=dot_product(b(i,:),b(iii,:))
1362 b(iii,:)=b(iii,:)-b(i,:)*maxp
1364 maxp=sqrt(dot_product(b(iii,:),b(iii,:)))
1365 b(iii,:)=b(iii,:)/maxp
1373 aa(ii, :)=aa(maxi,:)
1378 aa(ii,:)=aa(ii,:)/maxp
1381 maxp=dot_product(aa(ii,:),aa(i,:))
1383 aa(i,:)=aa(i,:)-aa(ii,:)*maxp
1401integer(spi),
intent(in ):: m,n
1402integer(spi),
dimension(n),
intent(in ):: ipiv
1403real(dp),
dimension(m,n),
intent(in ):: tt
1404real(dp),
dimension(m),
intent(in ):: v
1405real(dp),
dimension(m),
intent(out):: vv
1418 vv(i)=vv(i)-vv(j)*tt(i,j)
1444integer(spi),
intent(in ):: m,n
1445logical,
dimension(m,n),
intent(in ):: mask
1446real(dp),
dimension(m,n),
intent(in ):: a
1447real(dp),
dimension(m ),
intent(out):: d
1448real(dp),
dimension(m,n),
intent(out):: aa
1449real(dp),
dimension( n),
intent(out):: e
1450real(dp),
dimension(0:m+n,0:m+n):: g
1451real(dp),
dimension(0:m+n) :: h
1452integer(spi) :: i,j,k,nh
1457 if(mask(i,j))aa(i,j)=log(abs(a(i,j)))
1503 aa(i,j)=a(i,j)/(d(i)*e(j))
1519real(dp),
dimension(3,3),
intent(in ):: orth33
1520real(dp),
dimension(3),
intent(out):: ax3
1521real(dp),
dimension(3,3) :: plane
1522real(dp),
dimension(3) :: x,y,z
1524integer(spi),
dimension(1):: ii
1525integer(spi) :: i,j,k
1527do i=1,3; z(i)=dot_product(plane(:,i),plane(:,i));
enddo
1529k=ii(1); i=1+mod(k,3); j=1+mod(i,3)
1531s=
absv(ax3);
if(s==0)
return
1537ax3=ax3*atan2(dot_product(y,z),dot_product(x,z))
1549real(dp),
dimension(3),
intent(in ):: ax3
1550real(dp),
dimension(3,3),
intent(out):: orth33
1551real(dp),
dimension(3,3):: ax33
1563complex(dpc),
dimension(2,2),
intent(IN ):: cspin
1564real(dp),
dimension(0:3),
intent(OUT):: q
1565q(0)=real(cspin(1,1)); q(3)=aimag(cspin(1,1))
1566q(2)=real(cspin(2,1)); q(1)=aimag(cspin(2,1))
1576real(dp),
dimension(0:3),
intent(IN ):: q
1577complex(dpc),
dimension(2,2),
intent(OUT):: cspin
1578cspin(1,1)=cmplx( q(0), q(3))
1579cspin(2,1)=cmplx( q(2), q(1))
1580cspin(1,2)=cmplx(-q(2), q(1))
1581cspin(2,2)=cmplx( q(0),-q(3))
1592real(dp),
dimension(3,3),
intent(IN ):: rot
1593real(dp),
dimension(0:3),
intent(OUT):: q
1594real(dp),
dimension(3,3) :: t1,t2
1595real(dp),
dimension(3) :: u1,u2
1596real(dp) :: gamma,gammah,s,ss
1598integer(spi),
dimension(1):: ii
1601t1=rot;
do i=1,3; t1(i,i)=t1(i,i)-1; u1(i)=dot_product(t1(i,:),t1(i,:));
enddo
1602ii=maxloc(u1); j=ii(1); ss=u1(j)
1604 q=zero; q(0)=one;
return
1606t1(j,:)=t1(j,:)/sqrt(ss)
1613 t1(i,:)=t1(i,:)-dot_product(t1(1,:),t1(i,:))*t1(1,:)
1614 u1(i)=dot_product(t1(i,:),t1(i,:))
1622if(ss==zero)stop
'In rotov; invalid rot'
1623if(j/=2)t1(2,:)=t1(3,:)
1624t1(2,:)=t1(2,:)/sqrt(ss)
1626t1(3,1)=t1(1,2)*t1(2,3)-t1(1,3)*t1(2,2)
1627t1(3,2)=t1(1,3)*t1(2,1)-t1(1,1)*t1(2,3)
1628t1(3,3)=t1(1,1)*t1(2,2)-t1(1,2)*t1(2,1)
1630t2=matmul(t1,matmul(rot,transpose(t1)))
1632gamma=atan2(t2(2,1),t2(1,1)); gammah=gamma/two
1647real(dp),
dimension(0:3),
intent(IN ):: q
1648real(dp),
dimension(3,3),
intent(OUT):: rot
1649call setem(q(0),q(1),q(2),q(3),rot)
1659real(dp),
dimension(3),
intent(in ):: v
1660real(dp),
dimension(0:3),
intent(out):: q
1661real(dp),
dimension(3,3):: rot
1673real(dp),
dimension(0:3),
intent(in ):: q
1674real(dp),
dimension(3),
intent(out):: v
1675real(dp),
dimension(3,3):: rot
1691real(dp),
intent(IN ):: c,d,e,g
1692real(dp),
dimension(3,3),
intent(OUT):: r
1693real(dp):: cc,dd,ee,gg,de,dg,eg,dc,ec,gc
1694cc=c*c; dd=d*d; ee=e*e; gg=g*g
1695de=d*e; dg=d*g; eg=e*g
1696dc=d*c; ec=e*c; gc=g*c
1697r(1,1)=cc+dd-ee-gg; r(2,2)=cc-dd+ee-gg; r(3,3)=cc-dd-ee+gg
1698r(2,3)=2*(eg-dc); r(3,1)=2*(dg-ec); r(1,2)=2*(de-gc)
1699r(3,2)=2*(eg+dc); r(1,3)=2*(dg+ec); r(2,1)=2*(de+gc)
1710real(dp),
dimension(0:3),
intent(IN ):: a,b
1711real(dp),
dimension(0:3) :: c
1712c(0)=a(0)*b(0) -a(1)*b(1) -a(2)*b(2) -a(3)*b(3)
1713c(1)=a(0)*b(1) +a(1)*b(0) +a(2)*b(3) -a(3)*b(2)
1714c(2)=a(0)*b(2) +a(2)*b(0) +a(3)*b(1) -a(1)*b(3)
1715c(3)=a(0)*b(3) +a(3)*b(0) +a(1)*b(2) -a(2)*b(1)
1731integer(spi),
intent(IN ):: n
1732real(dp),
dimension(n,n),
intent(IN ):: a
1733real(dp),
dimension(n,n),
intent(OUT):: b
1734real(dp),
intent(OUT):: detb
1735integer(spi),
parameter :: L=5
1736real(dp),
dimension(n,n):: c,p
1739m=10+floor(log(
u1+maxval(abs(a)))/log(
u2))
1754detb=
u0;
do i=1,n; detb=detb+a(i,i); enddo; detb=exp(detb)
1769integer(spi),
intent(IN ):: n
1770real(dp),
dimension(n,n),
intent(IN ):: a
1771real(dp),
dimension(n,n),
intent(OUT):: b
1772real(dp),
dimension(n,n,(n*(n+1))/2),
intent(OUT):: bd
1773real(dp),
intent(OUT):: detb
1774real(dp),
dimension((n*(n+1))/2),
intent(OUT):: detbd
1775integer(spi),
parameter :: L=5
1776real(dp),
dimension(n,n) :: c,p
1777real(dp),
dimension(n,n,(n*(n+1))/2):: pd,cd
1779integer(spi) :: i,j,k,m,n1
1781m=10+floor(log(
u1+maxval(abs(a)))/log(
u2))
1797if(k/=n1)stop
'In expmatd; n1 is inconsistent with n'
1803 pd(:,:,k)=(matmul(cd(:,:,k),p)+matmul(c,pd(:,:,k)))/i
1811 bd(:,:,k)=2*bd(:,:,k)+matmul(bd(:,:,k),b)+matmul(b,bd(:,:,k))
1818detb=
u0;
do i=1,n; detb=detb+a(i,i); enddo; detb=exp(detb)
1819detbd=
u0;
do k=1,n; detbd(k)=detb;
enddo
1836integer(spi),
intent(IN ):: n
1837real(dp),
dimension(n,n),
intent(IN ):: a
1838real(dp),
dimension(n,n),
intent(OUT):: b
1839real(dp),
dimension(n,n,(n*(n+1))/2),
intent(OUT):: bd
1840real(dp),
dimension(n,n,(n*(n+1))/2,(n*(n+1))/2),
intent(OUT):: bdd
1841real(dp),
intent(OUT):: detb
1842real(dp),
dimension((n*(n+1))/2),
intent(OUT):: detbd
1843real(dp),
dimension((n*(n+1))/2,(n*(n+1))/2),
intent(OUT):: detbdd
1844integer(spi),
parameter :: L=5
1845real(dp),
dimension(n,n) :: c,p
1846real(dp),
dimension(n,n,(n*(n+1))/2) :: pd,cd
1847real(dp),
dimension(n,n,(n*(n+1))/2,(n*(n+1))/2):: pdd,cdd
1849integer(spi) :: i,j,k,ki,kj,m,n1
1851m=10+floor(log(
u1+maxval(abs(a)))/log(
u2))
1868if(k/=n1)stop
'In expmatd; n1 is inconsistent with n'
1877 pdd(:,:,ki,kj)=(matmul(cd(:,:,ki),pd(:,:,kj)) &
1878 + matmul(cd(:,:,kj),pd(:,:,ki)) &
1879 + matmul(c,pdd(:,:,ki,kj)))/i
1883 pd(:,:,k)=(matmul(cd(:,:,k),p)+matmul(c,pd(:,:,k)))/i
1893 bdd(:,:,ki,kj)=
u2*bdd(:,:,ki,kj) &
1894 +matmul(bdd(:,:,ki,kj),b) &
1895 +matmul(bd(:,:,ki),bd(:,:,kj)) &
1896 +matmul(bd(:,:,kj),bd(:,:,ki)) &
1897 +matmul(b,bdd(:,:,ki,kj))
1901 bd(:,:,k)=2*bd(:,:,k)+matmul(bd(:,:,k),b)+matmul(b,bd(:,:,k))
1908detb=
u0;
do i=1,n; detb=detb+a(i,i); enddo; detb=exp(detb)
1909detbd=
u0;
do k=1,n; detbd(k)=detb;
enddo
1910detbdd=
u0;
do ki=1,n;
do kj=1,n; detbdd(ki,kj)=detb; enddo;
enddo
1923integer(spi),
intent(IN ):: n
1924real(dp),
intent(IN ):: z
1925real(dp),
intent(OUT):: zn
1926integer(spi),
parameter:: ni=100
1927real(dp),
parameter :: eps0=1.e-16_dp
1928integer(spi) :: i,i2,n2
1941 t=t*z2/(i2*(i2+n2-1))
1943 if(abs(t)<eps)
return
1945print
'("In zntay; full complement of iterations used")'
1962integer(spi),
intent(IN ):: n
1963real(dp),
intent(IN ):: z
1964real(dp),
intent(OUT):: zn,znd,zndd,znddd
1966integer(spi):: i,i2p3
1971 call zntay(n+1,z,znd)
1972 call zntay(n+2,z,zndd)
1973 call zntay(n+3,z,znddd)
1979 znddd=(znd-
u3*zndd)/z2
1985 znddd=(znd-i2p3*zndd)/z2
1991 znddd=-(znd-
u3*zndd)/z2
1997 znddd=-(znd-i2p3*zndd)/z2
2029real(dp),
dimension(3),
intent(IN ):: v
2030complex(dpc),
intent(OUT):: z
2031logical,
intent(OUT):: infz
2034z=cmplx(v(1),v(2),dpc)
2039 infz=(rr==
u0);
if(infz)
return
2054complex(dpc),
intent(IN ):: z
2055logical,
intent(IN ):: infz
2056real(dp),
dimension(3),
intent(OUT):: v
2057real(dp),
parameter:: zero=0_dp,one=1_dp,two=2_dp
2058real(dp) :: r,q,rs,rsc,rsbi
2059if(infz)then; v=(/zero,zero,-one/); return;
endif
2060r=real(z); q=aimag(z); rs=r*r+q*q
2084complex(dpc),
intent(IN ):: z
2085logical,
intent(IN ):: infz
2086real(dp),
dimension(3),
intent(OUT):: v
2087complex(dpc),
dimension(3),
intent(OUT):: vd
2088real(dp),
parameter :: zero=0_dp,one=1_dp,two=2_dp,four=4_dp
2089real(dp) :: r,q,rs,rsc,rsbi,rsbis
2090real(dp),
dimension(3):: u1,u2
2092if(infz)then; v=(/zero,zero,-one/); return;
endif
2093r=real(z); q=aimag(z); rs=r*r+q*q
2100u1(1)=two*(one+q*q-r*r)*rsbis
2101u1(2)=-four*r*q*rsbis
2105 vd(i)=cmplx(u1(i),-u2(i),dpc)
2124real(dp),
dimension(3),
intent(IN ):: xc0,xc1,xc2
2125complex(dpc),
intent(OUT):: aa,bb,cc,dd
2126real(dp),
parameter:: zero=0_dp,one=1_dp
2127logical :: infz0,infz1,infz2
2128complex(dpc) :: z0,z1,z2,z02,z10,z21
2129call ctoz(xc0,z0,infz0)
2130call ctoz(xc1,z1,infz1)
2131call ctoz(xc2,z2,infz2)
2136if( (z0==z1.and.infz0.eqv.infz1).or.&
2137 (z1==z2.and.infz1.eqv.infz2).or.&
2138 (z2==z0.and.infz2.eqv.infz0)) &
2139 stop
'In setmobius; anchor points must be distinct'
2140if(infz2 .or. (.not.infz0 .and. abs(z0)<abs(z2)))
then
2147 aa=sqrt(-z21/(z02*z10))
2167 cc=sqrt(-z10/(z02*z21))
2207complex(dp),
intent(IN ):: z0,z1,z2
2208logical,
intent(IN ):: infz0,infz1,infz2
2209complex(dpc),
intent(OUT):: aa,bb,cc,dd
2210real(dp),
parameter:: zero=0_dp,one=1_dp
2211complex(dpc) :: z02,z10,z21
2215if( (z0==z1.and.infz0.eqv.infz1).or.&
2216 (z1==z2.and.infz1.eqv.infz2).or.&
2217 (z2==z0.and.infz2.eqv.infz0)) &
2218 stop
'In setmobius; anchor points must be distinct'
2219if(infz2 .or. (.not.infz0 .and. abs(z0)<abs(z2)))
then
2226 aa=sqrt(-z21/(z02*z10))
2246 cc=sqrt(-z10/(z02*z21))
2278complex(dpc),
intent(IN ):: aa,bb,cc,dd,z
2279logical,
intent(IN ):: infz
2280complex(dpc),
intent(OUT):: w
2281logical,
intent(OUT):: infw
2282real(dp),
parameter:: zero=0_dp
2283complex(dpc) :: top,bot
2294if(abs(bot)==zero)
then
2313complex(dpc),
intent(IN ):: aa,bb,cc,dd
2314real(dp),
dimension(3),
intent(IN ):: vz
2315real(dp),
dimension(3),
intent(OUT):: vw
2318call ctoz(vz, z,infz)
2319call zmobius(aa,bb,cc,dd, z,infz, w,infw)
2320call ztoc(w,infw, vw)
2337complex(dpc),
intent(IN ):: aa,bb,cc,dd,zz
2338logical,
intent(IN ):: infz
2339complex(dpc),
intent(OUT):: zw
2340logical,
intent(OUT):: infw
2341real(dp),
parameter:: one=1_dp
2342complex(dpc) :: aai,bbi,cci,ddi,d
2348call zmobius(aai,bbi,cci,ddi, zz,infz, zw,infw)
Some of the commonly used constants (pi etc) mainly for double-precision subroutines.
real(dp), parameter u1
one
real(dp), parameter o2
half
real(dp), parameter u3
three
logical, parameter t
for pain-relief in logical ops
real(dp), parameter u0
zero
real(dp), parameter u2
two
Standard integer, real, and complex single and double precision kinds.
integer, parameter dp
Double precision real kind.
integer, parameter dpc
Double precision real kind.
integer, parameter sp
Single precision real kind.
integer, parameter dpi
Double precision integer kind.
integer, parameter spi
Single precision integer kind.
Module for handy vector and matrix operations in Euclidean geometry.
real(sp) function, dimension(size(a), size(a)) diagn_s(a)
Return the diagonal matrix whose elements are the given vector.
real(dp) function, dimension(4) triple_cross_product_d(u, v, w)
Return the triple_cross_product for 4-vectors.
real(dp) function det_d(a)
Return the determinant of a double precision matrix.
real(dp) function absv_d(a)
Return the absolute magnitude of a double precision real vector.
integer(spi) function, dimension(size(a), size(a)) diagn_i(a)
Return the diagonal matrix whose elements are the given vector.
subroutine rowgram(m, n, a, ipiv, tt, b, rank)
Without changing (tall) rectangular input matrix a, perform pivoted gram- Schmidt operations to ortho...
subroutine graml_d(as, b, nrank, detsign, ldet)
A version of gram_d where the determinant information is returned in logarithmic form (to avoid overf...
subroutine expmatd(n, a, b, bd, detb, detbd)
Like expmat, but for the 1st derivatives also.
integer(dpi) function det_id(a)
Return the determinant of a double precision integer matrix.
real(sp) function, dimension(size(u)) orthogonalized_s(u, a)
Return the part of vector a that is orthogonal to unit vector u.
real(sp) function, dimension(3) cross_product_s(a, b)
Return the cross product of two single precision real 3-vectors.
real(sp) function, dimension(4) triple_cross_product_s(u, v, w)
Deliver the triple-cross-product, x, of the three 4-vectors, u, v, w, with the sign convention that o...
integer(spi) function, dimension(3, 3) identity3_i()
Return the 3-dimensional integer identity matrix.
subroutine zmobius(aa, bb, cc, dd, z, infz, w, infw)
Perform a complex Mobius transformation from (z,infz) to (w,infw) where the transformation coefficien...
real(dp) function, dimension(3) axial33_d(b)
Return the 3-vector corresponding to the given antisymmetric "axial vector" matrix,...
real(dp) function sarea_d(v1, v2, v3)
Compute the area of the spherical triangle, {v1,v2,v3}.
real(sp) function trace_s(b)
Return the trace of a given single precision real matrix.
real(sp) function, dimension(size(a)) normalized_s(a)
Return the normalized version of a single precision real vector.
real(dp) function triple_product_d(a, b, c)
Return the triple product of three double precision real 3-vectors.
integer(spi) function trace_i(b)
Return the trace of a given integer matrix.
subroutine expmatdd(n, a, b, bd, bdd, detb, detbd, detbdd)
Like expmat, but for the 1st and 2nd derivatives also.
integer(spi) function, dimension(size(a), size(b)) outer_product_i(a, b)
Return the outer product matrix of two integer vectors.
real(sp) function dqarea_s(rlat, drlata, drlona, drlatb, drlonb, drlatc, drlonc)
Compute the area of the spherical quadrilateral with a vertex at latitude rlat, and three other verti...
integer(spi) function, dimension(size(b, 1)) diagnn_i(b)
Return the vector whose elements are the diagonal ones of a given matrix.
real(dp) function huarea_d(sa, sb)
Spherical area of right-angle triangle whose orthogonal sides have orthographic projection dimensions...
real(dp) function, dimension(size(a), size(a)) diagn_d(a)
Return the diagonal matrix whose elements are the given vector.
subroutine dlltoxy_s(rlat, drlat, drlon, x2)
From a reference latitude, and increments of latitude and longitude, return the local cartesian 2-vec...
real(sp) function hav_s(t)
Haversine function in single precision.
real(dp) function, dimension(size(a), size(b)) outer_product_d(a, b)
Return the outer product matrix of two double precision real vectors.
real(sp) function absv_s(a)
Return the absolute magnitude of a single precision real vector.
real(dp) function, dimension(size(b, 1)) diagnn_d(b)
Return the vector whose elements are the diagonal ones of a given matrix.
real(dp) function, dimension(3, 3) axial3_d(a)
Return the axial "vector", as an antisymmetric matrix, corresponding to the given 3-vector assuming a...
real(sp) function det_s(a)
Return the determinant of a single precision matrix.
subroutine plaingram_d(b, nrank)
A "plain" (unpivoted) version of Gram-Schmidt, for square matrices only.
real(sp) function sarea_s(v1, v2, v3)
Compute the area of the spherical triangle, {v1,v2,v3}, measured in the right-handed sense,...
subroutine gram_s(as, b, nrank, det)
Apply a form of Gram-Schmidt orthogonalization process to return as many normalized orthogonal basis ...
integer(spi) function, dimension(n, n) identity_i(n)
Return the integer identity matrix for a given dimensionality.
subroutine gram_d(as, b, nrank, det)
Apply a form of Gram-Schmidt orthogonalization process to return as many normalized orthogonal basis ...
subroutine dlltoxy_d(rlat, drlat, drlon, x2)
From a reference latitude, and increments of latitude and longitude, return the local cartesian 2-vec...
subroutine cmobius(aa, bb, cc, dd, vz, vw)
Perform a complex Mobius transformation from cartesian vz to cartesian vw where the transformation co...
real(dp) function dqarea_d(rlat, drlata, drlona, drlatb, drlonb, drlatc, drlonc)
Compute the area of the spherical quadrilateral with a vertex at latitude rlat, and three other verti...
subroutine normalize_d(v)
Normalize the given double precision real vector.
subroutine normalize_s(v)
Normalize the given single precision real vector.
real(dp) function, dimension(3) cross_product_d(a, b)
Return the cross product of two double precision real 3-vectors.
subroutine ztocd(z, infz, v, vd)
The convention adopted for the complex derivative is that, for a complex infinitesimal map displaceme...
real(sp) function, dimension(size(b, 1)) diagnn_s(b)
Return the vector whose elements are the diagonal ones of a given matrix.
subroutine plaingram_s(b, nrank)
A "plain" (unpivoted) version of Gram-Schmidt, for square matrices only.
integer(spi) function det_i(a)
Return the determinant of a single precision integer matrix.
real(sp) function triple_product_s(a, b, c)
Return the triple product of three single precision real 3-vectors.
real(dp) function hav_d(t)
Haversine function in double precision.
subroutine zsetmobius(z0, infz0, z1, infz1, z2, infz2, aa, bb, cc, dd)
Find the Mobius transformation complex coefficients, aa,bb,cc,dd, with aa*dd-bb*cc=1,...
real(sp) function huarea_s(sa, sb)
Spherical area of right-angle triangle whose orthogonal sides have orthographic projection dimensions...
real(sp) function, dimension(3, 3) axial3_s(a)
Return the axial "vector", as an antisymmetric matrix, corresponding to the given 3-vector assuming a...
real(sp) function dtarea_s(rlat, drlata, drlona, drlatb, drlonb)
Compute the area of the spherical triangle with a vertex at latitude rlat, and two other vertices,...
real(dp) function, dimension(size(u)) orthogonalized_d(u, a)
Return the part of vector a that is orthogonal to unit vector u.
real(dp) function trace_d(b)
Return the trace of a given double precision real matrix.
real(sp) function, dimension(size(a), size(b)) outer_product_s(a, b)
Return the outer product matrix of two single precision real vectors.
real(dp) function dtarea_d(rlat, drlata, drlona, drlatb, drlonb)
Compute the area of the spherical triangle with a vertex at latitude rlat, and two other vertices,...
subroutine zmobiusi(aa, bb, cc, dd, zz, infz, zw, infw)
Perform the inverse of the mobius transformation with coefficients, {aa,bb,cc,dd}.
real(sp) function, dimension(3) axial33_s(b)
Return the 3-vector corresponding to the given antisymmetric "axial vector" matrix,...
real(dp) function, dimension(size(a)) normalized_d(a)
Return the normalized version of a double precision real vector.
subroutine cyclic(u1, u2, u3, d1, d2, d3)
Cyclically permute real vectors, u1, u2, u3, and scalars, d1, d2, d3.