grid_tools 1.14.0
Loading...
Searching...
No Matches
pmat4.f90
Go to the documentation of this file.
1
5
45module pmat4
46!============================================================================
47use pkind, only: spi,sp,dp,dpc
48implicit none
49private
54 axtoq,qtoax, &
59
60interface absv; module procedure absv_s,absv_d; end interface
61interface normalized;module procedure normalized_s,normalized_d;end interface
63 module procedure orthogonalized_s,orthogonalized_d; end interface
65 module procedure cross_product_s,cross_product_d, &
69 end interface
71 module procedure triple_product_s,triple_product_d; end interface
72interface det; module procedure det_s,det_d,det_i,det_id; end interface
73interface axial
74 module procedure axial3_s,axial3_d,axial33_s,axial33_d; end interface
75interface diag
77 end interface
78interface trace; module procedure trace_s,trace_d,trace_i; end interface
79interface identity; module procedure identity_i,identity3_i; end interface
80interface huarea; module procedure huarea_s,huarea_d; end interface
81interface sarea
83 end interface
84interface dlltoxy; module procedure dlltoxy_s,dlltoxy_d; end interface
85interface hav; module procedure hav_s, hav_d; end interface
86interface normalize;module procedure normalize_s,normalize_d; end interface
87interface gram
89 end interface
90interface rowops; module procedure rowops; end interface
91interface corral; module procedure corral; end interface
92interface rottoax; module procedure rottoax; end interface
93interface axtorot; module procedure axtorot; end interface
94interface spintoq; module procedure spintoq; end interface
95interface qtospin; module procedure qtospin; end interface
96interface rottoq; module procedure rottoq; end interface
97interface qtorot; module procedure qtorot; end interface
98interface axtoq; module procedure axtoq; end interface
99interface qtoax; module procedure qtoax; end interface
100interface setem; module procedure setem; end interface
101interface mulqq; module procedure mulqq; end interface
102interface expmat; module procedure expmat,expmatd,expmatdd; end interface
103interface zntay; module procedure zntay; end interface
104interface znfun; module procedure znfun; end interface
105interface ctoz; module procedure ctoz; end interface
106interface ztoc; module procedure ztoc,ztocd; end interface
107interface setmobius;module procedure setmobius,zsetmobius; end interface
108interface mobius; module procedure zmobius,cmobius; end interface
109interface mobiusi; module procedure zmobiusi; end interface
110
111contains
112
118function absv_s(a)result(s)! [absv]
119implicit none
120real(sp),dimension(:),intent(in):: a
121real(sp) :: s
122s=sqrt(dot_product(a,a))
123end function absv_s
124
130function absv_d(a)result(s)! [absv]
131implicit none
132real(dp),dimension(:),intent(in):: a
133real(dp) :: s
134s=sqrt(dot_product(a,a))
135end function absv_d
136
142function normalized_s(a)result(b)! [normalized]
143use pietc_s, only: u0
144implicit none
145real(sp),dimension(:),intent(IN):: a
146real(sp),dimension(size(a)) :: b
147real(sp) :: s
148s=absv_s(a); if(s==u0)then; b=u0;else;b=a/s;endif
149end function normalized_s
150
156function normalized_d(a)result(b)! [normalized]
157use pietc, only: u0
158implicit none
159real(dp),dimension(:),intent(IN):: a
160real(dp),dimension(size(a)) :: b
161real(dp) :: s
162s=absv_d(a); if(s==u0)then; b=u0;else;b=a/s;endif
163end function normalized_d
164
171function orthogonalized_s(u,a)result(b)! [orthogonalized]
172implicit none
173real(sp),dimension(:),intent(in):: u,a
174real(sp),dimension(size(u)) :: b
175real(sp) :: s
176! Note: this routine assumes u is already normalized
177s=dot_product(u,a); b=a-u*s
178end function orthogonalized_s
179
186function orthogonalized_d(u,a)result(b)! [orthogonalized]
187implicit none
188real(dp),dimension(:),intent(in):: u,a
189real(dp),dimension(size(u)) :: b
190real(dp) :: s
191! Note: this routine assumes u is already normalized
192s=dot_product(u,a); b=a-u*s
193end function orthogonalized_d
194
201function cross_product_s(a,b)result(c)! [cross_product]
202implicit none
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)
206end function cross_product_s
207
214function cross_product_d(a,b)result(c)! [cross_product]
215implicit none
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)
219end function cross_product_d
220
231function triple_cross_product_s(u,v,w)result(x)! [cross_product]
232implicit none
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)
243end function triple_cross_product_s
244
252function triple_cross_product_d(u,v,w)result(x)! [cross_product]
253implicit none
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)
264end function triple_cross_product_d
265
272function outer_product_s(a,b)result(c)! [outer_product]
273implicit none
274real(sp),dimension(:), intent(in ):: a
275real(sp),dimension(:), intent(in ):: b
276real(sp),DIMENSION(size(a),size(b)):: c
277integer(spi) :: nb,i
278nb=size(b)
279do i=1,nb; c(:,i)=a*b(i); enddo
280end function outer_product_s
281
288function outer_product_d(a,b)result(c)! [outer_product]
289implicit none
290real(dp),dimension(:), intent(in ):: a
291real(dp),dimension(:), intent(in ):: b
292real(dp),dimension(size(a),size(b)):: c
293integer(spi) :: nb,i
294nb=size(b)
295do i=1,nb; c(:,i)=a*b(i); enddo
296end function outer_product_d
297
304function outer_product_i(a,b)result(c)! [outer_product]
305implicit none
306integer(spi),dimension(:), intent(in ):: a
307integer(spi),dimension(:), intent(in ):: b
308integer(spi),dimension(size(a),size(b)):: c
309integer(spi) :: nb,i
310nb=size(b)
311do i=1,nb; c(:,i)=a*b(i); enddo
312end function outer_product_i
313
321function triple_product_s(a,b,c)result(tripleproduct)! [triple_product]
322implicit none
323real(sp),dimension(3),intent(IN ):: a,b,c
324real(sp) :: tripleproduct
325tripleproduct=dot_product( cross_product(a,b),c )
326end function triple_product_s
327
335function triple_product_d(a,b,c)result(tripleproduct)! [triple_product]
336implicit none
337real(dp),dimension(3),intent(IN ):: a,b,c
338real(dp) :: tripleproduct
339tripleproduct=dot_product( cross_product(a,b),c )
340end function triple_product_d
341
347function det_s(a)result(det)! [det]
348use pietc_s, only: u0
349implicit none
350real(sp),dimension(:,:),intent(IN ) :: a
351real(sp) :: det
352real(sp),dimension(size(a,1),size(a,1)):: b
353integer(spi) :: n,nrank
354n=size(a,1)
355if(n==3)then
356 det=triple_product(a(:,1),a(:,2),a(:,3))
357else
358 call gram(a,b,nrank,det)
359 if(nrank<n)det=u0
360endif
361end function det_s
362
368function det_d(a)result(det)! [det]
369use pietc, only: u0
370implicit none
371real(dp),dimension(:,:),intent(IN ) :: a
372real(dp) :: det
373real(dp),dimension(size(a,1),size(a,1)):: b
374integer(spi) :: n,nrank
375n=size(a,1)
376if(n==3)then
377 det=triple_product(a(:,1),a(:,2),a(:,3))
378else
379 call gram(a,b,nrank,det)
380 if(nrank<n)det=u0
381endif
382end function det_d
383
389function det_i(a)result(idet)! [det]
390implicit none
391integer(spi), dimension(:,:),intent(IN ):: a
392integer(spi) :: idet
393real(dp),dimension(size(a,1),size(a,2)):: b
394real(dp) :: bdet
395b=a; bdet=det(b); idet=nint(bdet)
396end function det_i
397
403function det_id(a)result(idet)! [det]
404use pkind, only: dp,dpi
405implicit none
406integer(dpi), dimension(:,:),intent(IN ):: a
407integer(dpi) :: idet
408real(dp),dimension(size(a,1),size(a,2)) :: b
409real(dp) :: bdet
410b=a; bdet=det(b); idet=nint(bdet)
411end function det_id
412
419function axial3_s(a)result(b)! [axial]
420use pietc_s, only: u0
421implicit none
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)
425end function axial3_s
426
433function axial3_d(a)result(b)! [axial]
434use pietc, only: u0
435implicit none
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)
439end function axial3_d
440
447function axial33_s(b)result(a)! [axial]
448use pietc_s, only: o2
449implicit none
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
453end function axial33_s
454
461function axial33_d(b)result(a)! [axial]
462use pietc, only: o2
463implicit none
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
467end function axial33_d
468
475function diagn_s(a)result(b)! [diag]
476use pietc, only: u0
477implicit none
478real(sp),dimension(:),intent(IN ) :: a
479real(sp),dimension(size(a),size(a)):: b
480integer(spi) :: n,i
481n=size(a)
482b=u0; do i=1,n; b(i,i)=a(i); enddo
483end function diagn_s
484
491function diagn_d(a)result(b)! [diag]
492use pietc, only: u0
493implicit none
494real(dp),dimension(:),intent(IN ) :: a
495real(dp),dimension(size(a),size(a)):: b
496integer(spi) :: n,i
497n=size(a)
498b=u0; do i=1,n; b(i,i)=a(i); enddo
499end function diagn_d
500
507function diagn_i(a)result(b)! [diag]
508implicit none
509integer(spi),dimension(:),intent(IN ) :: a
510integer(spi),dimension(size(a),size(a)):: b
511integer(spi) :: n,i
512n=size(a)
513b=0; do i=1,n; b(i,i)=a(i); enddo
514end function diagn_i
515
522function diagnn_s(b)result(a)! [diag]
523implicit none
524real(sp),dimension(:,:),intent(IN ):: b
525real(sp),dimension(size(b,1)) :: a
526integer(spi) :: n,i
527n=size(b,1)
528do i=1,n; a(i)=b(i,i); enddo
529end function diagnn_s
530
537function diagnn_d(b)result(a)! [diag]
538implicit none
539real(dp),dimension(:,:),intent(IN ):: b
540real(dp),dimension(size(b,1)) :: a
541integer(spi) :: n,i
542n=size(b,1)
543do i=1,n; a(i)=b(i,i); enddo
544end function diagnn_d
545
552function diagnn_i(b)result(a)! [diag]
553implicit none
554integer(spi),dimension(:,:),intent(IN ):: b
555integer(spi),dimension(size(b,1)) :: a
556integer(spi) :: n,i
557n=size(b,1)
558do i=1,n; a(i)=b(i,i); enddo
559end function diagnn_i
560
566function trace_s(b)result(s)! [trace]
567implicit none
568real(sp),dimension(:,:),intent(IN ):: b
569real(sp) :: s
570s=sum(diag(b))
571end function trace_s
572
578function trace_d(b)result(s)! [trace]
579implicit none
580real(dp),dimension(:,:),intent(IN ):: b
581real(dp) :: s
582s=sum(diag(b))
583end function trace_d
584
590function trace_i(b)result(s)! [trace]
591implicit none
592integer(spi),dimension(:,:),intent(IN ):: b
593integer(spi) :: s
594s=sum(diag(b))
595end function trace_i
596
602function identity_i(n)result(a)! [identity]
603implicit none
604integer(spi),intent(IN ) :: n
605integer(spi),dimension(n,n):: a
606integer(spi) :: i
607a=0; do i=1,n; a(i,i)=1; enddo
608end function identity_i
609
614function identity3_i()result(a)! [identity]
615implicit none
616integer(spi),dimension(3,3):: a
617integer(spi) :: i
618a=0; do i=1,3; a(i,i)=1; enddo
619end function identity3_i
620
629function huarea_s(sa,sb)result(area)! [huarea]
630implicit none
631real(sp),intent(IN ):: sa,sb
632real(sp) :: area
633real(sp) :: ca,cb
634ca=sqrt(1-sa**2)
635cb=sqrt(1-sb**2)
636area=asin(sa*sb/(1+ca*cb))
637end function huarea_s
638
647function huarea_d(sa,sb)result(area)! [huarea]
648implicit none
649real(dp),intent(IN ):: sa,sb
650real(dp) :: area
651real(dp) :: ca,cb
652ca=sqrt(1-sa**2)
653cb=sqrt(1-sb**2)
654area=asin(sa*sb/(1+ca*cb))
655end function huarea_d
656
667function sarea_s(v1,v2,v3)result(area)! [sarea]
668use pietc_s, only: zero=>u0
669implicit none
670real(sp),dimension(3),intent(IN ):: v1,v2,v3
671real(sp) :: area
672real(sp) :: s123,a1,a2,b,d1,d2,d3
673real(sp),dimension(3) :: u0,u1,u2,u3,x,y
674area=zero
675u1=normalized(v1); u2=normalized(v2); u3=normalized(v3)
676s123=triple_product(u1,u2,u3)
677if(s123==zero)return
678
679d1=dot_product(u3-u2,u3-u2)
680d2=dot_product(u1-u3,u1-u3)
681d3=dot_product(u2-u1,u2-u1)
682
683! Triangle that is not degenerate. Cyclically permute, so side 3 is longest:
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)
686y=normalized( cross_product(u1,u2) )
687b=dot_product(y,u3)
688u0=normalized( u3-y*b )
689x=cross_product(y,u0)
690a1=-dot_product(x,u1-u0); a2= dot_product(x,u2-u0)
691area=huarea(a1,b)+huarea(a2,b)
692
693contains
694
704subroutine cyclic(u1,u2,u3,d1,d2,d3)
705implicit none
706 real(sp),dimension(3),intent(INOUT):: u1,u2,u3
707 real(sp), intent(INOUT):: d1,d2,d3
708 real(sp),dimension(3) :: ut
709 real(sp) :: dt
710 dt=d1; d1=d2; d2=d3; d3=dt
711 ut=u1; u1=u2; u2=u3; u3=ut
712 end subroutine cyclic
713end function sarea_s
714
722function sarea_d(v1,v2,v3)result(area)! [sarea]
723use pietc, only: zero=>u0
724implicit none
725real(dp),dimension(3),intent(IN ):: v1,v2,v3
726real(dp) :: area
727real(dp) :: s123,a1,a2,b,d1,d2,d3
728real(dp),dimension(3) :: u0,u1,u2,u3,x,y
729area=zero
732if(s123==zero)return
733
734d1=dot_product(u3-u2,u3-u2)
735d2=dot_product(u1-u3,u1-u3)
736d3=dot_product(u2-u1,u2-u1)
737
738! Triangle that is not degenerate. Cyclically permute, so side 3 is longest:
739if(d3<d1 .or. d3<d2)call cyclic(u1,u2,u3,d1,d2,d3)
740if(d3<d1 .or. d3<d2)call cyclic(u1,u2,u3,d1,d2,d3)
742b=dot_product(y,u3)
743u0=normalized( u3-y*b )
744x=cross_product(y,u0)
745a1=-dot_product(x,u1-u0); a2= dot_product(x,u2-u0)
746area=huarea(a1,b)+huarea(a2,b)
747
748contains
749
759subroutine cyclic(u1,u2,u3,d1,d2,d3)
760implicit none
761 real(dp),dimension(3),intent(INOUT):: u1,u2,u3
762 real(dp), intent(INOUT):: d1,d2,d3
763 real(dp),dimension(3) :: ut
764 real(dp) :: dt
765 dt=d1; d1=d2; d2=d3; d3=dt
766 ut=u1; u1=u2; u2=u3; u3=ut
767 end subroutine cyclic
768end function sarea_d
769
785function dtarea_s(rlat,drlata,drlona,drlatb,drlonb) result(area)! [sarea]
786use pietc_s, only: u0,u1
787implicit none
788real(sp),intent(in ):: rlat,drlata,drlona,drlatb,drlonb
789real(sp) :: area
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
796cb=sqrt(u1-ssb)
797! Construct 2D normalized right-handed basis vectors with xb pointing to B:
798xb=x2b/sb
799yb=(/-xb(2),xb(1)/)
800xa=dot_product(xb,x2a)
801sa=dot_product(yb,x2a)
802ca=sqrt(u1-sa**2)
803sc=xa/ca
804cc=sqrt(u1-sc**2)
805sb=sb*cc-cb*sc
806area=huarea(-sa,sb)+huarea(sc,-sa)
807end function dtarea_s
808
824function dtarea_d(rlat,drlata,drlona,drlatb,drlonb) result(area)! [sarea]
825use pietc, only: u0,u1
826implicit none
827real(dp),intent(in ):: rlat,drlata,drlona,drlatb,drlonb
828real(dp) :: area
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
835cb=sqrt(u1-ssb)
836! Construct 2D normalized right-handed basis vectors with xb pointing to B:
837xb=x2b/sb
838yb=(/-xb(2),xb(1)/)
839xa=dot_product(xb,x2a)
840sa=dot_product(yb,x2a)
841ca=sqrt(u1-sa**2)
842sc=xa/ca
843cc=sqrt(u1-sc**2)
844sb=sb*cc-cb*sc
845area=huarea(-sa,sb)+huarea(sc,-sa)
846end function dtarea_d
847
867function dqarea_s &! [sarea]
868 (rlat,drlata,drlona,drlatb,drlonb,drlatc,drlonc) result(area)
869implicit none
870real(sp),intent(in ):: rlat,drlata,drlona,drlatb,drlonb,drlatc,drlonc
871real(sp) :: area
872area=sarea(rlat,drlata,drlona,drlatb,drlonb)&
873 -sarea(rlat,drlatc,drlonc,drlatb,drlonb)
874end function dqarea_s
875
876
896function dqarea_d &! [sarea]
897 (rlat,drlata,drlona,drlatb,drlonb,drlatc,drlonc) result(area)
898implicit none
899real(dp),intent(in ):: rlat,drlata,drlona,drlatb,drlonb,drlatc,drlonc
900real(dp) :: area
901area=sarea(rlat,drlata,drlona,drlatb,drlonb)&
902 -sarea(rlat,drlatc,drlonc,drlatb,drlonb)
903end function dqarea_d
904
915subroutine dlltoxy_s(rlat,drlat,drlon,x2)! [dlltoxy]
916use pietc_s, only: u2
917implicit none
918real(sp), intent(in ):: rlat,drlat,drlon
919real(sp),dimension(2),intent(out):: x2
920real(sp):: clata
921clata=cos(rlat+drlat)
922x2=(/clata*sin(drlon),sin(drlat)+u2*sin(rlat)*clata*hav(drlon)/)
923end subroutine dlltoxy_s
924
935subroutine dlltoxy_d(rlat,drlat,drlon,x2)! [dlltoxy]
936use pietc, only: u2
937implicit none
938real(dp), intent(in ):: rlat,drlat,drlon
939real(dp),dimension(2),intent(out):: x2
940real(dp):: clata
941clata=cos(rlat+drlat)
942x2=(/clata*sin(drlon),sin(drlat)+u2*sin(rlat)*clata*hav(drlon)/)
943end subroutine dlltoxy_d
944
950function hav_s(t) result(a)! [hav]
951use pietc_s, only: o2
952implicit none
953real(sp),intent(in ):: t
954real(sp) :: a
955a=(sin(t*o2))**2
956end function hav_s
957
963function hav_d(t) result(a)! [hav]
964use pietc, only: o2
965implicit none
966real(dp),intent(in ):: t
967real(dp) :: a
968a=(sin(t*o2))**2
969end function hav_d
970
975subroutine normalize_s(v)! [normalize]
976use pietc_s, only: u0,u1
977implicit none
978real(sp),dimension(:),intent(inout):: v
979real(sp) :: s
980s=absv(v); if(s==0)then; v=u0; v(1)=u1; else; v=v/s; endif
981end subroutine normalize_s
982
987subroutine normalize_d(v)! [normalize]
988use pietc, only: u0,u1
989implicit none
990real(dp),dimension(:),intent(inout):: v
991real(dp) :: s
992s=absv(v); if(s==u0)then; v=0; v(1)=u1; else; v=v/s; endif
993end subroutine normalize_d
994
1006subroutine gram_s(as,b,nrank,det)! [gram]
1007use pietc_s, only: u0,u1
1008implicit none
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
1020n=size(as,1)
1021m=size(as,2)
1022if(n/=size(b,1) .or. n/=size(b,2))stop 'In gram; incompatible dimensions'
1023a=as
1024b=identity(n)
1025det=u1
1026val=maxval(abs(a))
1027if(val==u0)then
1028 nrank=0
1029 return
1030endif
1031vcrit=val*crit
1032nrank=min(n,m)
1033do k=1,n
1034 if(k>nrank)exit
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)) )
1038 if(val<=vcrit)then
1039 nrank=k-1
1040 exit
1041 endif
1042 i=ii(1)
1043 j=ii(2)
1044 tv=b(:,j)
1045 b(:,j)=-b(:,k)
1046 b(:,k)=tv
1047 tv=a(:,i)
1048 a(:,i)=-a(:,k)
1049 a(:,k)=tv
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))
1053 s=sqrt(s)
1054 if(w(k)<u0)s=-s
1055 det=det*s
1056 b(:,k)=b(:,k)/s
1057 do l=k,n
1058 do j=l+1,n
1059 s=dot_product(b(:,l),b(:,j))
1060 b(:,j)=normalized( b(:,j)-b(:,l)*s )
1061 enddo
1062 enddo
1063enddo
1064end subroutine gram_s
1065
1077subroutine gram_d(as,b,nrank,det)! [gram]
1078use pietc, only: u0,u1
1079implicit none
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
1091n=size(as,1)
1092m=size(as,2)
1093if(n/=size(b,1) .or. n/=size(b,2))stop 'In gram; incompatible dimensions'
1094a=as
1095b=identity(n)
1096det=u1
1097val=maxval(abs(a))
1098if(val==u0)then
1099 nrank=0
1100 return
1101endif
1102vcrit=val*crit
1103nrank=min(n,m)
1104do k=1,n
1105 if(k>nrank)exit
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)) )
1109 if(val<=vcrit)then
1110 nrank=k-1
1111 exit
1112 endif
1113 i=ii(1)
1114 j=ii(2)
1115 tv=b(:,j)
1116 b(:,j)=-b(:,k)
1117 b(:,k)=tv
1118 tv=a(:,i)
1119 a(:,i)=-a(:,k)
1120 a(:,k)=tv
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))
1124 s=sqrt(s)
1125 if(w(k)<u0)s=-s
1126 det=det*s
1127 b(:,k)=b(:,k)/s
1128 do l=k,n
1129 do j=l+1,n
1130 s=dot_product(b(:,l),b(:,j))
1131 b(:,j)=normalized( b(:,j)-b(:,l)*s )
1132 enddo
1133 enddo
1134enddo
1135end subroutine gram_d
1136
1149subroutine graml_d(as,b,nrank,detsign,ldet)! [gram]
1150use pietc, only: u0
1151implicit none
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
1164detsign=1
1165n=size(as,1)
1166m=size(as,2)
1167if(n/=size(b,1) .or. n/=size(b,2))stop 'In gram; incompatible dimensions'
1168a=as
1169b=identity(n)
1170ldet=u0
1171val=maxval(abs(a))
1172if(val==u0)then
1173 nrank=0
1174 return
1175endif
1176vcrit=val*crit
1177nrank=min(n,m)
1178do k=1,n
1179 if(k>nrank)exit
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)) )
1183 if(val<=vcrit)then
1184 nrank=k-1
1185 exit
1186 endif
1187 i=ii(1)
1188 j=ii(2)
1189 tv=b(:,j)
1190 b(:,j)=-b(:,k)
1191 b(:,k)=tv
1192 tv=a(:,i)
1193 a(:,i)=-a(:,k)
1194 a(:,k)=tv
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))
1198 s=sqrt(s)
1199 if(w(k)<u0)s=-s
1200 if(s<0)then
1201 ldet=ldet+log(-s)
1202 detsign=-detsign
1203 elseif(s>u0)then
1204 ldet=ldet+log(s)
1205 else
1206 detsign=0
1207 endif
1208
1209 b(:,k)=b(:,k)/s
1210 do l=k,n
1211 do j=l+1,n
1212 s=dot_product(b(:,l),b(:,j))
1213 b(:,j)=normalized( b(:,j)-b(:,l)*s )
1214 enddo
1215 enddo
1216enddo
1217end subroutine graml_d
1218
1219
1226subroutine plaingram_s(b,nrank)! [gram]
1227use pietc_s, only: u0
1228implicit none
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'
1235val=maxval(abs(b))
1236nrank=0
1237if(val==0)then
1238 b=u0
1239 return
1240endif
1241vcrit=val*crit
1242do k=1,n
1243 val=sqrt(dot_product(b(:,k),b(:,k)))
1244 if(val<=vcrit)then
1245 b(:,k:n)=u0
1246 return
1247 endif
1248 b(:,k)=b(:,k)/val
1249 nrank=k
1250 do j=k+1,n
1251 b(:,j)=b(:,j)-b(:,k)*dot_product(b(:,k),b(:,j))
1252 enddo
1253enddo
1254end subroutine plaingram_s
1255
1262subroutine plaingram_d(b,nrank)! [gram]
1263use pietc, only: u0
1264implicit none
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'
1271val=maxval(abs(b))
1272nrank=0
1273if(val==u0)then
1274 b=u0
1275 return
1276endif
1277vcrit=val*crit
1278do k=1,n
1279 val=sqrt(dot_product(b(:,k),b(:,k)))
1280 if(val<=vcrit)then
1281 b(:,k:n)=u0
1282 return
1283 endif
1284 b(:,k)=b(:,k)/val
1285 nrank=k
1286 do j=k+1,n
1287 b(:,j)=b(:,j)-b(:,k)*dot_product(b(:,k),b(:,j))
1288 enddo
1289enddo
1290end subroutine plaingram_d
1291
1314subroutine rowgram(m,n,a,ipiv,tt,b,rank)! [gram]
1315use pietc, only: u0,u1
1316implicit none
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'
1331nepss=n*epss
1332rank=n
1333aa=a
1334tt=u0
1335do ii=1,n
1336! At this stage, all rows less than ii are already orthonormalized and are
1337! orthogonal to all rows at and beyond ii. Find the norms of these lower
1338! rows and pivot the largest of them into position ii:
1339 maxp=u0
1340 maxi=ii
1341 do i=ii,m
1342 p(i)=dot_product(aa(i,:),aa(i,:))
1343 if(p(i)>maxp)then
1344 maxp=p(i)
1345 maxi=i
1346 endif
1347 enddo
1348 if(maxp<nepss)then
1349 b=u0
1350 b(1:ii-1,:)=aa(1:ii-1,:)
1351! fill the remaining rows, ii:n, of b with remaining orthonormal rows at random
1352 do iii=ii,n
1353! find the column of b for which the maximum element is the smallest:
1354 do j=1,n
1355 rowv(j)=maxval(abs(b(1:iii-1,j)))
1356 enddo
1357 jloc=minloc(rowv)
1358 j=jloc(1)
1359 b(iii,j)=u1
1360 do i=1,iii-1
1361 maxp=dot_product(b(i,:),b(iii,:))
1362 b(iii,:)=b(iii,:)-b(i,:)*maxp
1363 enddo
1364 maxp=sqrt(dot_product(b(iii,:),b(iii,:)))
1365 b(iii,:)=b(iii,:)/maxp
1366 enddo
1367 rank=ii-1
1368 return
1369 endif
1370 ipiv(ii)=maxi
1371 if(maxi/=ii)then
1372 rowv =aa(ii, :)
1373 aa(ii, :)=aa(maxi,:)
1374 aa(maxi,:)=rowv
1375 endif
1376 maxp=sqrt(maxp)
1377 tt(ii,ii)=maxp
1378 aa(ii,:)=aa(ii,:)/maxp
1379! Adjust all rows below to make them orthogonal to new row ii
1380 do i=ii+1,m
1381 maxp=dot_product(aa(ii,:),aa(i,:))
1382 tt(i,ii)=maxp
1383 aa(i,:)=aa(i,:)-aa(ii,:)*maxp
1384 enddo
1385enddo
1386b=aa(1:n,:)
1387end subroutine rowgram
1388
1399subroutine rowops(m,n,ipiv,tt,v,vv)! [rowops]
1400implicit none
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
1406integer(spi):: i,j,k
1407real(dp) :: p
1408vv=v
1409do j=1,n
1410 k=ipiv(j)
1411 if(k/=j)then
1412 p=vv(j)
1413 vv(j)=vv(k)
1414 vv(k)=p
1415 endif
1416 vv(j)=vv(j)/tt(j,j)
1417 do i=j+1,m
1418 vv(i)=vv(i)-vv(j)*tt(i,j)
1419 enddo
1420enddo
1421end subroutine rowops
1422
1440subroutine corral(m,n,mask,a,d,aa,e)! [corral]
1441use pietc, only: u0,u1
1442use pmat, only: inv
1443implicit none
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
1453nh=1+m+n
1454aa=u0
1455do j=1,n
1456do i=1,m
1457 if(mask(i,j))aa(i,j)=log(abs(a(i,j)))
1458enddo
1459enddo
1460h=u0
1461g=u0
1462! Equations on row 0 enforcing Lagrange multiplier F-constraint:
1463do j=1,n
1464 k=m+j
1465 g(0,k)=u1
1466enddo
1467! Equations on rows 1:m minimizing row sums of quadratic terms:
1468do i=1,m
1469 do j=1,n
1470 k=m+j
1471 if(mask(i,j))then
1472 g(i,i)=g(i,i)-u1
1473 g(i,k)=-u1
1474 h(i)=h(i)-aa(i,j)
1475 endif
1476 enddo
1477enddo
1478! Equations on rows m+1:m+n minimizing col sums subject to constraint
1479do j=1,n
1480 k=m+j
1481 g(k,0)=u1
1482 do i=1,m
1483 if(mask(i,j))then
1484 g(k,k)=g(k,k)-u1
1485 g(k,i)=-u1
1486 h(k)=h(k)-aa(i,j)
1487 endif
1488 enddo
1489enddo
1490! Invert the normal equations:
1491call inv(g,h)
1492! Exponentiate the parts that become final scaling diagnonal matrices d and e:
1493do i=1,m
1494 d(i)=exp(h(i))
1495enddo
1496do j=1,n
1497 k=m+j
1498 e(j)=exp(h(k))
1499enddo
1500! Compute the rescaled matrix directly:
1501do j=1,n
1502do i=1,m
1503 aa(i,j)=a(i,j)/(d(i)*e(j))
1504enddo
1505enddo
1506end subroutine corral
1507
1517subroutine rottoax(orth33,ax3)! [rottoax]
1518implicit none
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
1523real(dp) :: s
1524integer(spi),dimension(1):: ii
1525integer(spi) :: i,j,k
1526plane=orth33-identity()! Columns must be coplanar vectors
1527do i=1,3; z(i)=dot_product(plane(:,i),plane(:,i)); enddo
1528ii=minloc(z)
1529k=ii(1); i=1+mod(k,3); j=1+mod(i,3)
1530ax3=cross_product(plane(:,i),plane(:,j))
1531s=absv(ax3); if(s==0)return
1532ax3=ax3/s ! <- temporarily a unit vector pointing along rotation axis
1533! Construct a unit 2D basis, x,y, in the plane of rotation
1534x=normalized(cross_product(ax3,plane(:,j)))
1535y=cross_product(ax3,x)
1536z=matmul(orth33,x)! Rotate x by the original rotation matrix
1537ax3=ax3*atan2(dot_product(y,z),dot_product(x,z))! multiply ax3 by the angle
1538end subroutine rottoax
1539
1547subroutine axtorot(ax3,orth33)! [axtorot]
1548implicit none
1549real(dp),dimension(3), intent(in ):: ax3
1550real(dp),dimension(3,3),intent(out):: orth33
1551real(dp),dimension(3,3):: ax33
1552real(dp) :: d
1553ax33=axial(ax3); call expmat(3,ax33,orth33,d)
1554end subroutine axtorot
1555
1561subroutine spintoq(cspin,q)! [spintoq]
1562implicit none
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))
1567end subroutine spintoq
1568
1574subroutine qtospin(q,cspin)! [qtospin]
1575implicit none
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))
1582end subroutine qtospin
1583
1589subroutine rottoq(rot,q)! [rottoq]
1590use pietc, only: zero=>u0,one=>u1,two=>u2
1591implicit none
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
1597integer(spi) :: i,j
1598integer(spi),dimension(1):: ii
1599! construct the orthogonal matrix, t1, whose third row is the rotation axis
1600! of rot:
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)
1603if(ss<1.e-16_dp)then
1604 q=zero; q(0)=one; return
1605endif
1606t1(j,:)=t1(j,:)/sqrt(ss)
1607if(j/=1)then
1608 u2 =t1(1,:)
1609 t1(1,:)=t1(j,:)
1610 t1(j,:)=u2
1611endif
1612do i=2,3
1613 t1(i,:)=t1(i,:)-dot_product(t1(1,:),t1(i,:))*t1(1,:)
1614 u1(i)=dot_product(t1(i,:),t1(i,:))
1615enddo
1616if(u1(3)>u1(2))then
1617 j=3
1618else
1619 j=2
1620endif
1621ss=u1(j)
1622if(ss==zero)stop 'In rotov; invalid rot'
1623if(j/=2)t1(2,:)=t1(3,:)
1624t1(2,:)=t1(2,:)/sqrt(ss)
1625! Form t1(3,:) as the cross product of t1(1,:) and t1(2,:)
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)
1629! Project rot into the frame whose axes are the rows of t1:
1630t2=matmul(t1,matmul(rot,transpose(t1)))
1631! Obtain the rotation angle, gamma, implied by rot, and gammah=gamma/2:
1632gamma=atan2(t2(2,1),t2(1,1)); gammah=gamma/two
1633! Hence deduce coefficients (in the form of a real 4-vector) of one of the two
1634! possible equivalent spinors:
1635s=sin(gammah)
1636q(0)=cos(gammah)
1637q(1:3)=t1(3,:)*s
1638end subroutine rottoq
1639
1645subroutine qtorot(q,rot)! [qtorot]
1646implicit none
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)
1650end subroutine qtorot
1651
1657subroutine axtoq(v,q)! [axtoq]
1658implicit none
1659real(dp),dimension(3), intent(in ):: v
1660real(dp),dimension(0:3),intent(out):: q
1661real(dp),dimension(3,3):: rot
1662call axtorot(v,rot)
1663call rottoq(rot,q)
1664end subroutine axtoq
1665
1671subroutine qtoax(q,v)! [qtoax]
1672implicit none
1673real(dp),dimension(0:3),intent(in ):: q
1674real(dp),dimension(3), intent(out):: v
1675real(dp),dimension(3,3):: rot
1676call qtorot(q,rot)
1677call rottoax(rot,v)
1678end subroutine qtoax
1679
1689subroutine setem(c,d,e,g,r)! [setem]
1690implicit none
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)
1700end subroutine setem
1701
1708function mulqq(a,b)result(c)! [mulqq]
1709implicit none
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)
1716end function mulqq
1717
1728subroutine expmat(n,a,b,detb)! [expmat]
1729use pietc, only: u0,u1,u2,o2
1730implicit none
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
1737real(dp) :: t
1738integer(spi) :: i,m
1739m=10+floor(log(u1+maxval(abs(a)))/log(u2))
1740t=o2**m
1741c=a*t
1742p=c
1743b=p
1744do i=2,l
1745 p=matmul(p,c)/i
1746 b=b+p
1747enddo
1748do i=1,m
1749 b=b*u2+matmul(b,b)
1750enddo
1751do i=1,n
1752 b(i,i)=b(i,i)+u1
1753enddo
1754detb=u0; do i=1,n; detb=detb+a(i,i); enddo; detb=exp(detb)
1755end subroutine expmat
1756
1766subroutine expmatd(n,a,b,bd,detb,detbd)! [expmat]
1767use pietc, only: u0,u1,u2,o2
1768implicit none
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
1778real(dp) :: t
1779integer(spi) :: i,j,k,m,n1
1780n1=(n*(n+1))*o2
1781m=10+floor(log(u1+maxval(abs(a)))/log(u2))
1782t=o2**m
1783c=a*t
1784p=c
1785pd=u0
1786do k=1,n
1787 pd(k,k,k)=t
1788enddo
1789k=n
1790do i=1,n-1
1791 do j=i+1,n
1792 k=k+1
1793 pd(i,j,k)=t
1794 pd(j,i,k)=t
1795 enddo
1796enddo
1797if(k/=n1)stop 'In expmatd; n1 is inconsistent with n'
1798cd=pd
1799b=p
1800bd=pd
1801do i=2,l
1802 do k=1,n1
1803 pd(:,:,k)=(matmul(cd(:,:,k),p)+matmul(c,pd(:,:,k)))/i
1804 enddo
1805 p=matmul(c,p)/i
1806 b=b+p
1807 bd=bd+pd
1808enddo
1809do i=1,m
1810 do k=1,n1
1811 bd(:,:,k)=2*bd(:,:,k)+matmul(bd(:,:,k),b)+matmul(b,bd(:,:,k))
1812 enddo
1813 b=b*u2+matmul(b,b)
1814enddo
1815do i=1,n
1816 b(i,i)=b(i,i)+u1
1817enddo
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
1820end subroutine expmatd
1821
1833subroutine expmatdd(n,a,b,bd,bdd,detb,detbd,detbdd)! [expmat]
1834use pietc, only: u0,u1,u2,o2
1835implicit none
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
1848real(dp) :: t
1849integer(spi) :: i,j,k,ki,kj,m,n1
1850n1=(n*(n+1))/2
1851m=10+floor(log(u1+maxval(abs(a)))/log(u2))
1852t=o2**m
1853c=a*t
1854p=c
1855pd=u0
1856pdd=u0
1857do k=1,n
1858 pd(k,k,k)=t
1859enddo
1860k=n
1861do i=1,n-1
1862 do j=i+1,n
1863 k=k+1
1864 pd(i,j,k)=t
1865 pd(j,i,k)=t
1866 enddo
1867enddo
1868if(k/=n1)stop 'In expmatd; n1 is inconsistent with n'
1869cd=pd
1870cdd=u0
1871b=p
1872bd=pd
1873bdd=u0
1874do i=2,l
1875 do ki=1,n1
1876 do kj=1,n1
1877 pdd(:,:,ki,kj)=(matmul(cd(:,:,ki),pd(:,:,kj)) &
1878 + matmul(cd(:,:,kj),pd(:,:,ki)) &
1879 + matmul(c,pdd(:,:,ki,kj)))/i
1880 enddo
1881 enddo
1882 do k=1,n1
1883 pd(:,:,k)=(matmul(cd(:,:,k),p)+matmul(c,pd(:,:,k)))/i
1884 enddo
1885 p=matmul(c,p)/i
1886 b=b+p
1887 bd=bd+pd
1888 bdd=bdd+pdd
1889enddo
1890do i=1,m
1891 do ki=1,n1
1892 do kj=1,n1
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))
1898 enddo
1899 enddo
1900 do k=1,n1
1901 bd(:,:,k)=2*bd(:,:,k)+matmul(bd(:,:,k),b)+matmul(b,bd(:,:,k))
1902 enddo
1903 b=b*u2+matmul(b,b)
1904enddo
1905do i=1,n
1906 b(i,i)=b(i,i)+u1
1907enddo
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
1911end subroutine expmatdd
1912
1920subroutine zntay(n,z,zn)! [zntay]
1921use pietc, only: u2
1922implicit none
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
1929real(dp) :: t,eps,z2
1930z2=z*u2
1931n2=n*2
1932t=1
1933do i=1,n
1934 t=t/(i*2-1)
1935enddo
1936eps=t*eps0
1937zn=t
1938t=t
1939do i=1,ni
1940 i2=i*2
1941 t=t*z2/(i2*(i2+n2-1))
1942 zn=zn+t
1943 if(abs(t)<eps)return
1944enddo
1945print'("In zntay; full complement of iterations used")'
1946end subroutine zntay
1947
1959subroutine znfun(n,z,zn,znd,zndd,znddd)! [znfun]
1960use pietc, only: u0,u2,u3
1961implicit none
1962integer(spi),intent(IN ):: n
1963real(dp), intent(IN ):: z
1964real(dp), intent(OUT):: zn,znd,zndd,znddd
1965real(dp) :: z2,rz2
1966integer(spi):: i,i2p3
1967z2=abs(z*u2)
1968rz2=sqrt(z2)
1969if(z2<u2)then
1970 call zntay(n ,z,zn)
1971 call zntay(n+1,z,znd)
1972 call zntay(n+2,z,zndd)
1973 call zntay(n+3,z,znddd)
1974else
1975 if(z>u0)then
1976 zn=cosh(rz2)
1977 znd=sinh(rz2)/rz2
1978 zndd=(zn-znd)/z2
1979 znddd=(znd-u3*zndd)/z2
1980 do i=1,n
1981 i2p3=i*2+3
1982 zn=znd
1983 znd=zndd
1984 zndd=znddd
1985 znddd=(znd-i2p3*zndd)/z2
1986 enddo
1987 else
1988 zn=cos(rz2)
1989 znd=sin(rz2)/rz2
1990 zndd=-(zn-znd)/z2
1991 znddd=-(znd-u3*zndd)/z2
1992 do i=1,n
1993 i2p3=i*2+3
1994 zn=znd
1995 znd=zndd
1996 zndd=znddd
1997 znddd=-(znd-i2p3*zndd)/z2
1998 enddo
1999 endif
2000endif
2001end subroutine znfun
2002
2018
2026subroutine ctoz(v, z,infz)! [ctoz]
2027use pietc, only: u0,u1
2028implicit none
2029real(dp),dimension(3),intent(IN ):: v
2030complex(dpc), intent(OUT):: z
2031logical, intent(OUT):: infz
2032real(dp) :: rr,zzpi
2033infz=.false.
2034z=cmplx(v(1),v(2),dpc)
2035if(v(3)>u0)then
2036 zzpi=u1/(u1+v(3))
2037else
2038 rr=v(1)**2+v(2)**2
2039 infz=(rr==u0); if(infz)return ! <- The point is mapped to infinity (90S)
2040 zzpi=(u1-v(3))/rr
2041endif
2042z=z*zzpi
2043end subroutine ctoz
2044
2052subroutine ztoc(z,infz, v)! [ztoc]
2053implicit none
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
2061rsc=one-rs
2062rsbi=one/(one+rs)
2063v(1)=two*rsbi*r
2064v(2)=two*rsbi*q
2065v(3)=rsc*rsbi
2066end subroutine ztoc
2067
2082subroutine ztocd(z,infz, v,vd)! [ztoc]
2083implicit none
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
2091integer(spi) :: i
2092if(infz)then; v=(/zero,zero,-one/); return; endif
2093r=real(z); q=aimag(z); rs=r*r+q*q
2094rsc=one-rs
2095rsbi=one/(one+rs)
2096rsbis=rsbi**2
2097v(1)=two*rsbi*r
2098v(2)=two*rsbi*q
2099v(3)=rsc*rsbi
2100u1(1)=two*(one+q*q-r*r)*rsbis
2101u1(2)=-four*r*q*rsbis
2102u1(3)=-four*r*rsbis
2103u2=cross_product(v,u1)
2104do i=1,3
2105 vd(i)=cmplx(u1(i),-u2(i),dpc)
2106enddo
2107end subroutine ztocd
2108
2122subroutine setmobius(xc0,xc1,xc2, aa,bb,cc,dd)! [setmobius]
2123implicit none
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)
2132z21=z2-z1
2133z02=z0-z2
2134z10=z1-z0
2135
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
2141! z0 is finite and smaller than z2:
2142 if(infz1)then
2143 aa=one/sqrt(z02) ! <- z1 is infinite
2144 elseif(infz2)then
2145 aa=one/sqrt(z10) ! <- z2 is infinite
2146 else
2147 aa=sqrt(-z21/(z02*z10)) ! <- all zs are finite
2148 endif
2149 bb=-z0*aa
2150 if(infz1)then
2151 cc=aa ! <- z1 is infinite
2152 dd=-z2*aa !
2153 elseif(infz2)then
2154 cc=zero ! <- z2 is infinite
2155 dd=z10*aa !
2156 else
2157 cc=-(z10/z21)*aa ! <- all zs are finite
2158 dd= z2*(z10/z21)*aa !
2159 endif
2160else
2161! z2 is finite and smaller than z0:
2162 if(infz0)then
2163 cc=one/sqrt(z21) ! <- z0 is inifinite
2164 elseif(infz1)then
2165 cc=one/sqrt(z02) ! <- z1 is infinite
2166 else
2167 cc=sqrt(-z10/(z02*z21)) ! <- all zs are finite
2168 endif
2169 dd=-z2*cc
2170 if(infz0)then
2171 aa=zero ! <- z0 is inifinite
2172 bb=-z21*cc !
2173 elseif(infz1)then
2174 aa=cc ! <- z1 is infinite
2175 bb=-z0*cc !
2176 else
2177 aa=(-z21/z10)*cc ! <- all zs are finite
2178 bb=z0*(z21/z10)*cc !
2179 endif
2180endif
2181end subroutine setmobius
2182
2205subroutine zsetmobius(z0,infz0, z1,infz1, z2,infz2, aa,bb,cc,dd)
2206implicit none
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
2212z21=z2-z1
2213z02=z0-z2
2214z10=z1-z0
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
2220! z0 is finite and smaller than z2:
2221 if(infz1)then
2222 aa=one/sqrt(z02) ! <- z1 is infinite
2223 elseif(infz2)then
2224 aa=one/sqrt(z10) ! <- z2 is infinite
2225 else
2226 aa=sqrt(-z21/(z02*z10)) ! <- all zs are finite
2227 endif
2228 bb=-z0*aa
2229 if(infz1)then
2230 cc=aa ! <- z1 is infinite
2231 dd=-z2*aa !
2232 elseif(infz2)then
2233 cc=zero ! <- z2 is infinite
2234 dd=z10*aa !
2235 else
2236 cc=-(z10/z21)*aa ! <- all zs are finite
2237 dd= z2*(z10/z21)*aa !
2238 endif
2239else
2240! z2 is finite and smaller than z0:
2241 if(infz0)then
2242 cc=one/sqrt(z21) ! <- z0 is inifinite
2243 elseif(infz1)then
2244 cc=one/sqrt(z02) ! <- z1 is infinite
2245 else
2246 cc=sqrt(-z10/(z02*z21)) ! <- all zs are finite
2247 endif
2248 dd=-z2*cc
2249 if(infz0)then
2250 aa=zero ! <- z0 is inifinite
2251 bb=-z21*cc !
2252 elseif(infz1)then
2253 aa=cc ! <- z1 is infinite
2254 bb=-z0*cc !
2255 else
2256 aa=(-z21/z10)*cc ! <- all zs are finite
2257 bb=z0*(z21/z10)*cc !
2258 endif
2259endif
2260end subroutine zsetmobius
2261
2276subroutine zmobius(aa,bb,cc,dd, z,infz, w,infw)! [mobius]
2277implicit none
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
2284w=0
2285infw=.false.
2286if(infz)then
2287 top=aa
2288 bot=cc
2289else
2290 top=aa*z+bb
2291 bot=cc*z+dd
2292endif
2293
2294if(abs(bot)==zero)then
2295 infw=.true.
2296else
2297 w=top/bot
2298endif
2299end subroutine zmobius
2300
2311subroutine cmobius(aa,bb,cc,dd, vz,vw)! [mobius]
2312implicit none
2313complex(dpc), intent(IN ):: aa,bb,cc,dd
2314real(dp),dimension(3),intent(IN ):: vz
2315real(dp),dimension(3),intent(OUT):: vw
2316complex(dpc):: z,w
2317logical :: infz,infw
2318call ctoz(vz, z,infz)
2319call zmobius(aa,bb,cc,dd, z,infz, w,infw)
2320call ztoc(w,infw, vw)
2321end subroutine cmobius
2322
2335subroutine zmobiusi(aa,bb,cc,dd, zz,infz, zw,infw) ! [mobiusi]
2336implicit none
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
2343d=one/(aa*dd-bb*cc)
2344aai=dd*d
2345ddi=aa*d
2346bbi=-bb*d
2347cci=-cc*d
2348call zmobius(aai,bbi,cci,ddi, zz,infz, zw,infw)
2349end subroutine zmobiusi
2350
2351
2352end module pmat4
Some of the commonly used constants (pi etc) mainly for double-precision subroutines.
Definition pietc.f90:14
real(dp), parameter u1
one
Definition pietc.f90:20
real(dp), parameter o2
half
Definition pietc.f90:32
real(dp), parameter u3
three
Definition pietc.f90:24
logical, parameter t
for pain-relief in logical ops
Definition pietc.f90:17
real(dp), parameter u0
zero
Definition pietc.f90:19
real(dp), parameter u2
two
Definition pietc.f90:22
Standard integer, real, and complex single and double precision kinds.
Definition pkind.f90:7
integer, parameter dp
Double precision real kind.
Definition pkind.f90:11
integer, parameter dpc
Double precision real kind.
Definition pkind.f90:13
integer, parameter sp
Single precision real kind.
Definition pkind.f90:10
integer, parameter dpi
Double precision integer kind.
Definition pkind.f90:9
integer, parameter spi
Single precision integer kind.
Definition pkind.f90:8
Module for handy vector and matrix operations in Euclidean geometry.
Definition pmat4.f90:45
real(sp) function, dimension(size(a), size(a)) diagn_s(a)
Return the diagonal matrix whose elements are the given vector.
Definition pmat4.f90:476
real(dp) function, dimension(4) triple_cross_product_d(u, v, w)
Return the triple_cross_product for 4-vectors.
Definition pmat4.f90:253
real(dp) function det_d(a)
Return the determinant of a double precision matrix.
Definition pmat4.f90:369
real(dp) function absv_d(a)
Return the absolute magnitude of a double precision real vector.
Definition pmat4.f90:131
integer(spi) function, dimension(size(a), size(a)) diagn_i(a)
Return the diagonal matrix whose elements are the given vector.
Definition pmat4.f90:508
subroutine rowgram(m, n, a, ipiv, tt, b, rank)
Without changing (tall) rectangular input matrix a, perform pivoted gram- Schmidt operations to ortho...
Definition pmat4.f90:1315
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...
Definition pmat4.f90:1150
subroutine expmatd(n, a, b, bd, detb, detbd)
Like expmat, but for the 1st derivatives also.
Definition pmat4.f90:1767
integer(dpi) function det_id(a)
Return the determinant of a double precision integer matrix.
Definition pmat4.f90:404
real(sp) function, dimension(size(u)) orthogonalized_s(u, a)
Return the part of vector a that is orthogonal to unit vector u.
Definition pmat4.f90:172
real(sp) function, dimension(3) cross_product_s(a, b)
Return the cross product of two single precision real 3-vectors.
Definition pmat4.f90:202
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...
Definition pmat4.f90:232
integer(spi) function, dimension(3, 3) identity3_i()
Return the 3-dimensional integer identity matrix.
Definition pmat4.f90:615
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...
Definition pmat4.f90:2277
real(dp) function, dimension(3) axial33_d(b)
Return the 3-vector corresponding to the given antisymmetric "axial vector" matrix,...
Definition pmat4.f90:462
real(dp) function sarea_d(v1, v2, v3)
Compute the area of the spherical triangle, {v1,v2,v3}.
Definition pmat4.f90:723
real(sp) function trace_s(b)
Return the trace of a given single precision real matrix.
Definition pmat4.f90:567
real(sp) function, dimension(size(a)) normalized_s(a)
Return the normalized version of a single precision real vector.
Definition pmat4.f90:143
real(dp) function triple_product_d(a, b, c)
Return the triple product of three double precision real 3-vectors.
Definition pmat4.f90:336
integer(spi) function trace_i(b)
Return the trace of a given integer matrix.
Definition pmat4.f90:591
subroutine expmatdd(n, a, b, bd, bdd, detb, detbd, detbdd)
Like expmat, but for the 1st and 2nd derivatives also.
Definition pmat4.f90:1834
integer(spi) function, dimension(size(a), size(b)) outer_product_i(a, b)
Return the outer product matrix of two integer vectors.
Definition pmat4.f90:305
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...
Definition pmat4.f90:869
integer(spi) function, dimension(size(b, 1)) diagnn_i(b)
Return the vector whose elements are the diagonal ones of a given matrix.
Definition pmat4.f90:553
real(dp) function huarea_d(sa, sb)
Spherical area of right-angle triangle whose orthogonal sides have orthographic projection dimensions...
Definition pmat4.f90:648
real(dp) function, dimension(size(a), size(a)) diagn_d(a)
Return the diagonal matrix whose elements are the given vector.
Definition pmat4.f90:492
subroutine dlltoxy_s(rlat, drlat, drlon, x2)
From a reference latitude, and increments of latitude and longitude, return the local cartesian 2-vec...
Definition pmat4.f90:916
real(sp) function hav_s(t)
Haversine function in single precision.
Definition pmat4.f90:951
real(dp) function, dimension(size(a), size(b)) outer_product_d(a, b)
Return the outer product matrix of two double precision real vectors.
Definition pmat4.f90:289
real(sp) function absv_s(a)
Return the absolute magnitude of a single precision real vector.
Definition pmat4.f90:119
real(dp) function, dimension(size(b, 1)) diagnn_d(b)
Return the vector whose elements are the diagonal ones of a given matrix.
Definition pmat4.f90:538
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...
Definition pmat4.f90:434
real(sp) function det_s(a)
Return the determinant of a single precision matrix.
Definition pmat4.f90:348
subroutine plaingram_d(b, nrank)
A "plain" (unpivoted) version of Gram-Schmidt, for square matrices only.
Definition pmat4.f90:1263
real(sp) function sarea_s(v1, v2, v3)
Compute the area of the spherical triangle, {v1,v2,v3}, measured in the right-handed sense,...
Definition pmat4.f90:668
subroutine gram_s(as, b, nrank, det)
Apply a form of Gram-Schmidt orthogonalization process to return as many normalized orthogonal basis ...
Definition pmat4.f90:1007
integer(spi) function, dimension(n, n) identity_i(n)
Return the integer identity matrix for a given dimensionality.
Definition pmat4.f90:603
subroutine gram_d(as, b, nrank, det)
Apply a form of Gram-Schmidt orthogonalization process to return as many normalized orthogonal basis ...
Definition pmat4.f90:1078
subroutine dlltoxy_d(rlat, drlat, drlon, x2)
From a reference latitude, and increments of latitude and longitude, return the local cartesian 2-vec...
Definition pmat4.f90:936
subroutine cmobius(aa, bb, cc, dd, vz, vw)
Perform a complex Mobius transformation from cartesian vz to cartesian vw where the transformation co...
Definition pmat4.f90:2312
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...
Definition pmat4.f90:898
subroutine normalize_d(v)
Normalize the given double precision real vector.
Definition pmat4.f90:988
subroutine normalize_s(v)
Normalize the given single precision real vector.
Definition pmat4.f90:976
real(dp) function, dimension(3) cross_product_d(a, b)
Return the cross product of two double precision real 3-vectors.
Definition pmat4.f90:215
subroutine ztocd(z, infz, v, vd)
The convention adopted for the complex derivative is that, for a complex infinitesimal map displaceme...
Definition pmat4.f90:2083
real(sp) function, dimension(size(b, 1)) diagnn_s(b)
Return the vector whose elements are the diagonal ones of a given matrix.
Definition pmat4.f90:523
subroutine plaingram_s(b, nrank)
A "plain" (unpivoted) version of Gram-Schmidt, for square matrices only.
Definition pmat4.f90:1227
integer(spi) function det_i(a)
Return the determinant of a single precision integer matrix.
Definition pmat4.f90:390
real(sp) function triple_product_s(a, b, c)
Return the triple product of three single precision real 3-vectors.
Definition pmat4.f90:322
real(dp) function hav_d(t)
Haversine function in double precision.
Definition pmat4.f90:964
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,...
Definition pmat4.f90:2206
real(sp) function huarea_s(sa, sb)
Spherical area of right-angle triangle whose orthogonal sides have orthographic projection dimensions...
Definition pmat4.f90:630
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...
Definition pmat4.f90:420
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,...
Definition pmat4.f90:786
real(dp) function, dimension(size(u)) orthogonalized_d(u, a)
Return the part of vector a that is orthogonal to unit vector u.
Definition pmat4.f90:187
real(dp) function trace_d(b)
Return the trace of a given double precision real matrix.
Definition pmat4.f90:579
real(sp) function, dimension(size(a), size(b)) outer_product_s(a, b)
Return the outer product matrix of two single precision real vectors.
Definition pmat4.f90:273
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,...
Definition pmat4.f90:825
subroutine zmobiusi(aa, bb, cc, dd, zz, infz, zw, infw)
Perform the inverse of the mobius transformation with coefficients, {aa,bb,cc,dd}.
Definition pmat4.f90:2336
real(sp) function, dimension(3) axial33_s(b)
Return the 3-vector corresponding to the given antisymmetric "axial vector" matrix,...
Definition pmat4.f90:448
real(dp) function, dimension(size(a)) normalized_d(a)
Return the normalized version of a double precision real vector.
Definition pmat4.f90:157
subroutine cyclic(u1, u2, u3, d1, d2, d3)
Cyclically permute real vectors, u1, u2, u3, and scalars, d1, d2, d3.
Definition pmat4.f90:705