grid_tools 1.14.0
Loading...
Searching...
No Matches
pmat5.f90
Go to the documentation of this file.
1
4
7module cstgeo ! Constants for orientation and stretching of map
8use pkind, only: sp
9implicit none
10real(sp),dimension(3,3):: rotm
11real(sp) :: sc
12real(sp) :: sci
13end module cstgeo
14
17module dcstgeo ! Constants for orientation and stretching of map
18use pkind, only: dp
19implicit none
20real(dp),dimension(3,3):: rotm
21real(dp) :: sc
22real(dp) :: sci
23end module dcstgeo
24
27module pmat5
28use pkind, only: spi,sp,dp
29implicit none
30private
34interface ininmap; module procedure sininmap,dininmap; end interface
35interface inivmap; module procedure sinivmap,dinivmap; end interface
36interface ctogr; module procedure sctogr, dctogr; end interface
37interface grtoc
38 module procedure sgrtoc,dgrtoc, sgrtocd,dgrtocd, sgrtocdd,dgrtocdd
39 end interface
40interface ctog; module procedure sctog, dctog; end interface
41interface gtoc
42 module procedure sgtoc,dgtoc, sgtocd,dgtocd, sgtocdd,dgtocdd; end interface
43interface gtoframe
44 module procedure sgtoframev,gtoframev,sgtoframem,gtoframem; end interface
45interface paraframe; module procedure sparaframe,paraframe; end interface
46interface frametwist;module procedure sframetwist,frametwist; end interface
47interface ctoc_schm
48 module procedure sctoc,dctoc, sctocd,dctocd, sctocdd,dctocdd; end interface
49interface plrot; module procedure plrot, dplrot; end interface
50interface plroti; module procedure plroti,dplroti; end interface
51interface plctoc; module procedure plctoc; end interface
52contains
53
63subroutine sininmap(alon0,alat0,rot3)! [ininmap]
64use pietc_s, only: u0,dtor
65implicit none
66real(sp), intent(IN ):: alon0,alat0
67real(sp),dimension(3,3),intent(OUT):: rot3
68real(sp) :: blon0,blat0,clon0,clat0,slon0,slat0
69blon0=dtor*alon0; clon0=cos(blon0); slon0=sin(blon0)
70blat0=dtor*alat0; clat0=cos(blat0); slat0=sin(blat0)
71rot3(1,1)=slat0*clon0; rot3(1,2)=slat0*slon0; rot3(1,3)=-clat0
72rot3(2,1)=-slon0; rot3(2,2)=clon0; rot3(2,3)=u0
73rot3(3,1)=clat0*clon0; rot3(3,2)=clat0*slon0; rot3(3,3)=slat0
74end subroutine sininmap
75
85subroutine dininmap(alon0,alat0,rot3)! [ininmap]
86use pietc, only: u0,dtor
87implicit none
88real(dp), intent(IN ):: alon0,alat0
89real(dp),dimension(3,3),intent(OUT):: rot3
90real(dp) :: blon0,blat0,clon0,clat0,slon0,slat0
91blon0=dtor*alon0; clon0=cos(blon0); slon0=sin(blon0)
92blat0=dtor*alat0; clat0=cos(blat0); slat0=sin(blat0)
93rot3(1,1)=slat0*clon0; rot3(1,2)=slat0*slon0; rot3(1,3)=-clat0
94rot3(2,1)=-slon0; rot3(2,2)=clon0; rot3(2,3)=u0
95rot3(3,1)=clat0*clon0; rot3(3,2)=clat0*slon0; rot3(3,3)=slat0
96end subroutine dininmap
97
107subroutine sinivmap(alon0,alat0,rot3)! [inivmap]
108use pietc_s, only: u0,dtor
109implicit none
110real(sp), intent(IN ):: alon0,alat0
111real(sp),dimension(3,3),intent(OUT):: rot3
112real(sp) :: blon0,blat0,clon0,clat0,slon0,slat0
113blon0=dtor*alon0
114blat0=dtor*alat0
115clon0=cos(blon0)
116slon0=sin(blon0)
117clat0=cos(blat0)
118slat0=sin(blat0)
119rot3(1,1)=-slon0
120rot3(1,2)=clon0
121rot3(1,3)=u0
122rot3(2,1)=-slat0*clon0
123rot3(2,2)=-slat0*slon0
124rot3(2,3)=clat0
125rot3(3,1)=clat0*clon0
126rot3(3,2)=clat0*slon0
127rot3(3,3)=slat0
128end subroutine sinivmap
129
139subroutine dinivmap(alon0,alat0,rot3)! [inivmap]
140use pietc, only: u0,dtor
141implicit none
142real(dp), intent(IN ):: alon0,alat0
143real(dp),dimension(3,3),intent(OUT):: rot3
144real(dp) :: blon0,blat0,clon0,clat0,slon0,slat0
145blon0=dtor*alon0
146blat0=dtor*alat0
147clon0=cos(blon0)
148slon0=sin(blon0)
149clat0=cos(blat0)
150slat0=sin(blat0)
151rot3(1,1)=-slon0
152rot3(1,2)=clon0
153rot3(1,3)=u0
154rot3(2,1)=-slat0*clon0
155rot3(2,2)=-slat0*slon0
156rot3(2,3)=clat0
157rot3(3,1)=clat0*clon0
158rot3(3,2)=clat0*slon0
159rot3(3,3)=slat0
160end subroutine dinivmap
161
172subroutine sctogr(xe,rlat,rlon)! [ctogr]
173use pietc_s, only: u0
174implicit none
175real(sp),dimension(3),intent(IN ):: xe
176real(sp), intent(OUT):: rlat,rlon
177real(sp) :: r
178r=sqrt(xe(1)**2+xe(2)**2)
179rlat=atan2(xe(3),r)
180if(r==u0)then
181 rlon=u0
182else
183 rlon=atan2(xe(2),xe(1))
184endif
185end subroutine sctogr
186
197subroutine dctogr(xe,rlat,rlon)! [ctogr]
198use pietc, only: u0
199implicit none
200real(dp),dimension(3),intent(IN ):: xe
201real(dp), intent(OUT):: rlat,rlon
202real(dp) :: r
203r=sqrt(xe(1)**2+xe(2)**2)
204rlat=atan2(xe(3),r)
205if(r==u0)then
206 rlon=u0
207else
208 rlon=atan2(xe(2),xe(1))
209endif
210end subroutine dctogr
211
219subroutine sgrtoc(rlat,rlon,xe)! [grtoc]
220implicit none
221real(sp), intent(IN ):: rlat,rlon
222real(sp),dimension(3),intent(OUT):: xe
223real(sp) :: sla,cla,slo,clo
224sla=sin(rlat); cla=cos(rlat)
225slo=sin(rlon); clo=cos(rlon)
226xe(1)=cla*clo; xe(2)=cla*slo; xe(3)=sla
227end subroutine sgrtoc
228
236subroutine dgrtoc(rlat,rlon,xe)! [grtoc]
237implicit none
238real(dp), intent(IN ):: rlat,rlon
239real(dp),dimension(3),intent(OUT):: xe
240real(dp) :: sla,cla,slo,clo
241sla=sin(rlat); cla=cos(rlat)
242slo=sin(rlon); clo=cos(rlon)
243xe(1)=cla*clo; xe(2)=cla*slo; xe(3)=sla
244end subroutine dgrtoc
245
256subroutine sgrtocd(rlat,rlon,xe,dxedlat,dxedlon)! [grtoc]
257implicit none
258real(sp), intent(IN ):: rlat,rlon
259real(sp),dimension(3),intent(OUT):: xe,dxedlat,dxedlon
260real(dp) :: rlat_d,rlon_d
261real(dp),dimension(3):: xe_d,dxedlat_d,dxedlon_d
262rlat_d=rlat; rlon_d=rlon
263call dgrtocd(rlat_d,rlon_d,xe_d,dxedlat_d,dxedlon_d)
264xe =xe_d
265dxedlat=dxedlat_d
266dxedlon=dxedlon_d
267end subroutine sgrtocd
268
279subroutine dgrtocd(rlat,rlon,xe,dxedlat,dxedlon)! [grtoc]
280use pietc, only: u0
281implicit none
282real(dp), intent(IN ):: rlat,rlon
283real(dp),dimension(3),intent(OUT):: xe,dxedlat,dxedlon
284real(dp) :: sla,cla,slo,clo
285sla=sin(rlat); cla=cos(rlat)
286slo=sin(rlon); clo=cos(rlon)
287xe(1)=cla*clo; xe(2)=cla*slo; xe(3)=sla
288dxedlat(1)=-sla*clo; dxedlat(2)=-sla*slo; dxedlat(3)=cla
289dxedlon(1)=-cla*slo; dxedlon(2)= cla*clo; dxedlon(3)=u0
290end subroutine dgrtocd
291
305subroutine sgrtocdd(rlat,rlon,xe,dxedlat,dxedlon, &! [grtoc]
306 ddxedlatdlat,ddxedlatdlon,ddxedlondlon)
307implicit none
308real(sp), intent(IN ):: rlat,rlon
309real(sp),dimension(3),intent(OUT):: xe,dxedlat,dxedlon, &
310 ddxedlatdlat,ddxedlatdlon,ddxedlondlon
311real(dp) :: rlat_d,rlon_d
312real(dp),dimension(3):: xe_d,dxedlat_d,dxedlon_d, &
313 ddxedlatdlat_d,ddxedlatdlon_d,ddxedlondlon_d
314rlat_d=rlat; rlon_d=rlon
315call dgtocdd(rlat_d,rlon_d,xe_d,dxedlat_d,dxedlon_d, &
316 ddxedlatdlat_d,ddxedlatdlon_d,ddxedlondlon_d)
317xe =xe_d
318dxedlat =dxedlat_d
319dxedlon =dxedlon_d
320ddxedlatdlat=ddxedlatdlat_d
321ddxedlatdlon=ddxedlatdlon_d
322ddxedlondlon=ddxedlondlon_d
323end subroutine sgrtocdd
324
338subroutine dgrtocdd(rlat,rlon,xe,dxedlat,dxedlon, &! [grtoc]
339 ddxedlatdlat,ddxedlatdlon,ddxedlondlon)
340use pietc, only: u0
341implicit none
342real(dp), intent(IN ):: rlat,rlon
343real(dp),dimension(3),intent(OUT):: xe,dxedlat,dxedlon, &
344 ddxedlatdlat,ddxedlatdlon,ddxedlondlon
345real(dp) :: sla,cla,slo,clo
346sla=sin(rlat); cla=cos(rlat)
347slo=sin(rlon); clo=cos(rlon)
348xe(1)=cla*clo; xe(2)=cla*slo; xe(3)=sla
349dxedlat(1)=-sla*clo; dxedlat(2)=-sla*slo; dxedlat(3)=cla
350dxedlon(1)=-cla*slo; dxedlon(2)= cla*clo; dxedlon(3)=u0
351ddxedlatdlat(1)=-cla*clo
352ddxedlatdlat(2)=-cla*slo
353ddxedlatdlat(3)=-sla
354ddxedlatdlon(1)= sla*slo
355ddxedlatdlon(2)=-sla*clo
356ddxedlatdlon(3)= u0
357ddxedlondlon(1)=-cla*clo
358ddxedlondlon(2)=-cla*slo
359ddxedlondlon(3)= u0
360end subroutine dgrtocdd
361
372subroutine sctog(xe,dlat,dlon)! [ctog]
373use pietc_s, only: u0,rtod
374implicit none
375real(sp),dimension(3),intent(IN ):: xe
376real(sp), intent(OUT):: dlat,dlon
377real(sp) :: r
378r=sqrt(xe(1)**2+xe(2)**2)
379dlat=atan2(xe(3),r)*rtod
380if(r==u0)then
381 dlon=u0
382else
383 dlon=atan2(xe(2),xe(1))*rtod
384endif
385end subroutine sctog
386
397subroutine dctog(xe,dlat,dlon)! [ctog]
398use pietc, only: u0,rtod
399implicit none
400real(dp),dimension(3),intent(IN ):: xe
401real(dp), intent(OUT):: dlat,dlon
402real(dp) :: r
403r=sqrt(xe(1)**2+xe(2)**2)
404dlat=atan2(xe(3),r)*rtod
405if(r==u0)then
406 dlon=u0
407else
408 dlon=atan2(xe(2),xe(1))*rtod
409endif
410end subroutine dctog
411
422subroutine sgtoc(dlat,dlon,xe)! [gtoc]
423use pietc_s, only: dtor
424implicit none
425real(sp), intent(IN ):: dlat,dlon
426real(sp),dimension(3),intent(OUT):: xe
427real(sp) :: rlat,rlon,sla,cla,slo,clo
428rlat=dtor*dlat; rlon=dtor*dlon
429sla=sin(rlat); cla=cos(rlat)
430slo=sin(rlon); clo=cos(rlon)
431xe(1)=cla*clo; xe(2)=cla*slo; xe(3)=sla
432end subroutine sgtoc
433
444subroutine dgtoc(dlat,dlon,xe)! [gtoc]
445use pietc, only: dtor
446implicit none
447real(dp), intent(IN ):: dlat,dlon
448real(dp),dimension(3),intent(OUT):: xe
449real(dp) :: rlat,rlon,sla,cla,slo,clo
450rlat=dtor*dlat; rlon=dtor*dlon
451sla=sin(rlat); cla=cos(rlat)
452slo=sin(rlon); clo=cos(rlon)
453xe(1)=cla*clo; xe(2)=cla*slo; xe(3)=sla
454end subroutine dgtoc
455
469subroutine sgtocd(dlat,dlon,xe,dxedlat,dxedlon)! [gtoc]
470implicit none
471real(sp), intent(IN ):: dlat,dlon
472real(sp),dimension(3),intent(OUT):: xe,dxedlat,dxedlon
473real(dp) :: dlat_d,dlon_d
474real(dp),dimension(3):: xe_d,dxedlat_d,dxedlon_d
475dlat_d=dlat; dlon_d=dlon
476call dgtocd(dlat_d,dlon_d,xe_d,dxedlat_d,dxedlon_d)
477xe =xe_d
478dxedlat=dxedlat_d
479dxedlon=dxedlon_d
480end subroutine sgtocd
481
495subroutine dgtocd(dlat,dlon,xe,dxedlat,dxedlon)! [gtoc]
496use pietc, only: u0,dtor
497implicit none
498real(dp), intent(IN ):: dlat,dlon
499real(dp),dimension(3),intent(OUT):: xe,dxedlat,dxedlon
500real(dp) :: rlat,rlon,sla,cla,slo,clo
501rlat=dtor*dlat; rlon=dtor*dlon
502sla=sin(rlat); cla=cos(rlat)
503slo=sin(rlon); clo=cos(rlon)
504xe(1)=cla*clo; xe(2)=cla*slo; xe(3)=sla
505dxedlat(1)=-sla*clo; dxedlat(2)=-sla*slo; dxedlat(3)=cla; dxedlat=dxedlat*dtor
506dxedlon(1)=-cla*slo; dxedlon(2)= cla*clo; dxedlon(3)=u0 ; dxedlon=dxedlon*dtor
507end subroutine dgtocd
508
526subroutine sgtocdd(dlat,dlon,xe,dxedlat,dxedlon, &
527 ddxedlatdlat,ddxedlatdlon,ddxedlondlon)! [gtoc]
528implicit none
529real(sp), intent(IN ):: dlat,dlon
530real(sp),dimension(3),intent(OUT):: xe,dxedlat,dxedlon, &
531 ddxedlatdlat,ddxedlatdlon,ddxedlondlon
532real(dp) :: dlat_d,dlon_d
533real(dp),dimension(3):: xe_d,dxedlat_d,dxedlon_d, &
534 ddxedlatdlat_d,ddxedlatdlon_d,ddxedlondlon_d
535dlat_d=dlat; dlon_d=dlon
536call dgtocdd(dlat_d,dlon_d,xe_d,dxedlat_d,dxedlon_d, &
537 ddxedlatdlat_d,ddxedlatdlon_d,ddxedlondlon_d)
538xe =xe_d
539dxedlat =dxedlat_d
540dxedlon =dxedlon_d
541ddxedlatdlat=ddxedlatdlat_d
542ddxedlatdlon=ddxedlatdlon_d
543ddxedlondlon=ddxedlondlon_d
544end subroutine sgtocdd
545
563subroutine dgtocdd(dlat,dlon,xe,dxedlat,dxedlon, &
564 ddxedlatdlat,ddxedlatdlon,ddxedlondlon)! [gtoc]
565use pietc, only: u0,dtor
566implicit none
567real(dp), intent(IN ):: dlat,dlon
568real(dp),dimension(3),intent(OUT):: xe,dxedlat,dxedlon, &
569 ddxedlatdlat,ddxedlatdlon,ddxedlondlon
570real(dp) :: rlat,rlon,sla,cla,slo,clo
571rlat=dtor*dlat; rlon=dtor*dlon
572sla=sin(rlat); cla=cos(rlat)
573slo=sin(rlon); clo=cos(rlon)
574xe(1)=cla*clo; xe(2)=cla*slo; xe(3)=sla
575dxedlat(1)=-sla*clo; dxedlat(2)=-sla*slo; dxedlat(3)=cla; dxedlat=dxedlat*dtor
576dxedlon(1)=-cla*slo; dxedlon(2)= cla*clo; dxedlon(3)=u0 ; dxedlon=dxedlon*dtor
577ddxedlatdlat(1)=-cla*clo
578ddxedlatdlat(2)=-cla*slo
579ddxedlatdlat(3)=-sla
580ddxedlatdlon(1)= sla*slo
581ddxedlatdlon(2)=-sla*clo
582ddxedlatdlon(3)= u0
583ddxedlondlon(1)=-cla*clo
584ddxedlondlon(2)=-cla*slo
585ddxedlondlon(3)= u0
586ddxedlatdlat=ddxedlatdlat*dtor**2
587ddxedlatdlon=ddxedlatdlon*dtor**2
588ddxedlondlon=ddxedlondlon*dtor**2
589end subroutine dgtocdd
590
599subroutine sgtoframem(splat,splon,sorth)! [gtoframe]
600implicit none
601real(sp), intent(in ):: splat,splon
602real(sp),dimension(3,3),intent(out):: sorth
603real(dp):: plat,plon
604real(dp),dimension(3,3):: orth
605plat=splat; plon=splon; call gtoframem(plat,plon,orth); sorth=orth
606end subroutine sgtoframem
607
616subroutine gtoframem(plat,plon,orth)! [gtoframe]
617implicit none
618real(dp), intent(in ):: plat,plon
619real(dp),dimension(3,3),intent(out):: orth
620real(dp),dimension(3):: xp,yp,zp
621call gtoframev(plat,plon, xp,yp,zp)
622orth(:,1)=xp; orth(:,2)=yp; orth(:,3)=zp
623end subroutine gtoframem
624
636subroutine sgtoframev(splat,splon,sxp,syp,szp)! [gtoframe]
637!==============================================================================
638implicit none
639real(sp), intent(in ):: splat,splon
640real(sp),dimension(3),intent(out):: sxp,syp,szp
641!------------------------------------------------------------------------------
642real(dp) :: plat,plon
643real(dp),dimension(3):: xp,yp,zp
644!==============================================================================
645plat=splat; plon=splon
646call gtoframev(plat,plon, xp,yp,zp)
647sxp=xp; syp=yp; szp=zp
648end subroutine sgtoframev
649
661subroutine gtoframev(plat,plon, xp,yp,zp)! [gtoframe]
662use pietc, only: u0,u1
663implicit none
664real(dp), intent(in ):: plat,plon
665real(dp),dimension(3),intent(out):: xp,yp,zp
666real(dp),dimension(3):: dzpdlat,dzpdlon
667if(plat==90)then ! is this the north pole?
668 xp=(/ u0,u1, u0/) ! Assume the limiting case lat-->90 along the 0-meridian
669 yp=(/-u1,u0, u0/) !
670 zp=(/ u0,u0, u1/)
671elseif(plat==-90)then
672 xp=(/ u0,u1, u0/) ! Assume the limiting case lat-->90 along the 0-meridian
673 yp=(/ u1,u0, u0/) !
674 zp=(/ u0,u0,-u1/)
675else
676 call gtoc(plat,plon,zp,dzpdlat,dzpdlon)
677 xp=dzpdlon/sqrt(dot_product(dzpdlon,dzpdlon))
678 yp=dzpdlat/sqrt(dot_product(dzpdlat,dzpdlat))
679endif
680end subroutine gtoframev
681
695subroutine sparaframe(sxp,syp,szp, sxv,syv,szv)! [paraframe]
696implicit none
697real(sp),dimension(3),intent(in ):: sxp,syp,szp, szv
698real(sp),dimension(3),intent(out):: sxv,syv
699real(dp),dimension(3):: xp,yp,zp, xv,yv,zv
700xp=sxp; yp=syp; zp=szp
701call paraframe(xp,yp,zp, xv,yv,zv)
702sxv=xv; syv=yv
703end subroutine sparaframe
704
718subroutine paraframe(xp,yp,zp, xv,yv,zv)! [paraframe]
720implicit none
721real(dp),dimension(3),intent(in ):: xp,yp,zp, zv
722real(dp),dimension(3),intent(out):: xv,yv
723real(dp) :: xpofv,ypofv,theta,ctheta,stheta
724real(dp),dimension(3):: xq,yq
725xpofv=dot_product(xp,zv)
726ypofv=dot_product(yp,zv)
727theta=atan2(ypofv,xpofv); ctheta=cos(theta); stheta=sin(theta)
728xq=zv-zp; xq=xq-zv*dot_product(xq,zv); xq=normalized(xq)
729yq=cross_product(zv,xq)
730xv=xq*ctheta-yq*stheta
731yv=xq*stheta+yq*ctheta
732end subroutine paraframe
733
750subroutine sframetwist(sxp,syp,szp, sxv,syv,szv, stwist)! [frametwist]
751implicit none
752real(sp),dimension(3),intent(in ):: sxp,syp,szp, sxv,syv,szv
753real(sp), intent(out):: stwist
754real(dp),dimension(3):: xp,yp,zp, xv,yv,zv
755real(dp) :: twist
756xp=sxp;yp=syp; zp=szp; xv=sxv; yv=syv; zv=szv
757call frametwist(xp,yp,zp, xv,yv,zv, twist)
758stwist=twist
759end subroutine sframetwist
760
777subroutine frametwist(xp,yp,zp, xv,yv,zv, twist)! [frametwist]
778implicit none
779real(dp),dimension(3),intent(in ):: xp,yp,zp, xv,yv,zv
780real(dp), intent(out):: twist
781real(dp),dimension(3):: xxv,yyv
782real(dp) :: c,s
783call paraframe(xp,yp,zp, xxv,yyv,zv)
784c=dot_product(xv,xxv); s=dot_product(xv,yyv)
785twist=atan2(s,c)
786end subroutine frametwist
787
794subroutine sctoc(s,xc1,xc2)! [ctoc_schm]
795use pietc_s, only: u1,u2
796implicit none
797real(sp), intent(IN ):: s
798real(sp),dimension(3),intent(INOUT):: xc1,xc2
799real(sp) :: x,y,z,a,b,d,e,ab2,aa,bb,di,aapbb,aambb
800x=xc1(1); y=xc1(2); z=xc1(3)
801a=s+u1
802b=s-u1
803ab2=a*b*u2
804aa=a*a
805bb=b*b
806aapbb=aa+bb
807aambb=aa-bb
808d=aapbb-ab2*z
809e=aapbb*z-ab2
810di=u1/d
811xc2(1)=(aambb*x)*di
812xc2(2)=(aambb*y)*di
813xc2(3)=e*di
814end subroutine sctoc
815
824subroutine sctocd(s,xc1,xc2,dxc2)! [ctoc_schm]
825use pietc_s, only: u0,u1,u2
826implicit none
827real(sp),intent(IN) :: s
828real(sp),dimension(3), intent(INOUT):: xc1,xc2
829real(sp),dimension(3,3),intent( OUT):: dxc2
830real(sp) :: x,y,z,a,b,d,e, &
831 ab2,aa,bb,di,ddi,aapbb,aambb
832x=xc1(1); y=xc1(2); z=xc1(3)
833a=s+u1
834b=s-u1
835ab2=a*b*u2
836aa=a*a
837bb=b*b
838aapbb=aa+bb
839aambb=aa-bb
840d=aapbb-ab2*z
841e=aapbb*z-ab2
842di=u1/d
843xc2(1)=(aambb*x)*di
844xc2(2)=(aambb*y)*di
845xc2(3)=e*di
846ddi=di*di
847
848dxc2=u0
849dxc2(1,1)=aambb*di
850dxc2(1,3)=ab2*aambb*x*ddi
851dxc2(2,2)=aambb*di
852dxc2(2,3)=ab2*aambb*y*ddi
853dxc2(3,3)=aapbb*di +ab2*e*ddi
854end subroutine sctocd
855
865subroutine sctocdd(s,xc1,xc2,dxc2,ddxc2)! [ctoc_schm]
866use pietc_s, only: u0,u1,u2
867implicit none
868real(sp), intent(IN ):: s
869real(sp),dimension(3), intent(INOUT):: xc1,xc2
870real(sp),dimension(3,3), intent( OUT):: dxc2
871real(sp),dimension(3,3,3),intent( OUT):: ddxc2
872real(sp) :: x,y,z,a,b,d,e, &
873 ab2,aa,bb,di,ddi,dddi, &
874 aapbb,aambb
875x=xc1(1); y=xc1(2); z=xc1(3)
876a=s+u1
877b=s-u1
878ab2=a*b*u2
879aa=a*a
880bb=b*b
881aapbb=aa+bb
882aambb=aa-bb
883d=aapbb-ab2*z
884e=aapbb*z-ab2
885di=u1/d
886xc2(1)=(aambb*x)*di
887xc2(2)=(aambb*y)*di
888xc2(3)=e*di
889ddi=di*di
890dddi=ddi*di
891
892dxc2=u0
893dxc2(1,1)=aambb*di
894dxc2(1,3)=ab2*aambb*x*ddi
895dxc2(2,2)=aambb*di
896dxc2(2,3)=ab2*aambb*y*ddi
897dxc2(3,3)=aapbb*di +ab2*e*ddi
898
899ddxc2=u0
900ddxc2(1,1,3)=ab2*aambb*ddi
901ddxc2(1,3,1)=ddxc2(1,1,3)
902ddxc2(1,3,3)=u2*ab2**2*aambb*x*dddi
903ddxc2(2,2,3)=ab2*aambb*ddi
904ddxc2(2,3,2)=ddxc2(2,2,3)
905ddxc2(2,3,3)=u2*ab2**2*aambb*y*dddi
906ddxc2(3,3,3)=u2*ab2*(aapbb*ddi+ab2*e*dddi)
907end subroutine sctocdd
908
915subroutine dctoc(s,xc1,xc2)! [ctoc_schm]
916use pietc, only: u1,u2
917implicit none
918real(dp), intent(IN ):: s
919real(dp),dimension(3),intent(INOUT):: xc1,xc2
920real(dp) :: x,y,z,a,b,d,e, &
921 ab2,aa,bb,di,aapbb,aambb
922x=xc1(1); y=xc1(2); z=xc1(3)
923a=s+u1
924b=s-u1
925ab2=a*b*u2
926aa=a*a
927bb=b*b
928aapbb=aa+bb
929aambb=aa-bb
930d=aapbb-ab2*z
931e=aapbb*z-ab2
932di=u1/d
933xc2(1)=(aambb*x)*di
934xc2(2)=(aambb*y)*di
935xc2(3)=e*di
936end subroutine dctoc
937
946subroutine dctocd(s,xc1,xc2,dxc2)! [ctoc_schm]
947use pietc, only: u0,u1,u2
948implicit none
949real(dp), intent(IN ):: s
950real(dp),dimension(3), intent(INOUT):: xc1,xc2
951real(dp),dimension(3,3),intent( OUT):: dxc2
952real(dp) :: x,y,z,a,b,d,e, &
953 ab2,aa,bb,di,ddi,aapbb,aambb
954x=xc1(1); y=xc1(2); z=xc1(3)
955a=s+u1
956b=s-u1
957ab2=a*b*u2
958aa=a*a
959bb=b*b
960aapbb=aa+bb
961aambb=aa-bb
962d=aapbb-ab2*z
963e=aapbb*z-ab2
964di=u1/d
965xc2(1)=(aambb*x)*di
966xc2(2)=(aambb*y)*di
967xc2(3)=e*di
968ddi=di*di
969
970dxc2=u0
971dxc2(1,1)=aambb*di
972dxc2(1,3)=ab2*aambb*x*ddi
973dxc2(2,2)=aambb*di
974dxc2(2,3)=ab2*aambb*y*ddi
975dxc2(3,3)=aapbb*di +ab2*e*ddi
976end subroutine dctocd
977
987subroutine dctocdd(s,xc1,xc2,dxc2,ddxc2)! [ctoc_schm]
988use pietc, only: u0,u1,u2
989implicit none
990real(dp),intent(IN) :: s
991real(dp),dimension(3), intent(INOUT):: xc1,xc2
992real(dp),dimension(3,3), intent(OUT ):: dxc2
993real(dp),dimension(3,3,3),intent(OUT ):: ddxc2
994real(dp) :: x,y,z,a,b,d,e, &
995 ab2,aa,bb,di,ddi,dddi, &
996 aapbb,aambb
997x=xc1(1); y=xc1(2); z=xc1(3)
998a=s+u1
999b=s-u1
1000ab2=a*b*u2
1001aa=a*a
1002bb=b*b
1003aapbb=aa+bb
1004aambb=aa-bb
1005d=aapbb-ab2*z
1006e=aapbb*z-ab2
1007di=u1/d
1008xc2(1)=(aambb*x)*di
1009xc2(2)=(aambb*y)*di
1010xc2(3)=e*di
1011ddi=di*di
1012dddi=ddi*di
1013
1014dxc2=u0
1015dxc2(1,1)=aambb*di
1016dxc2(1,3)=ab2*aambb*x*ddi
1017dxc2(2,2)=aambb*di
1018dxc2(2,3)=ab2*aambb*y*ddi
1019dxc2(3,3)=aapbb*di +ab2*e*ddi
1020
1021ddxc2=u0
1022ddxc2(1,1,3)=ab2*aambb*ddi
1023ddxc2(1,3,1)=ddxc2(1,1,3)
1024ddxc2(1,3,3)=u2*ab2**2*aambb*x*dddi
1025ddxc2(2,2,3)=ab2*aambb*ddi
1026ddxc2(2,3,2)=ddxc2(2,2,3)
1027ddxc2(2,3,3)=u2*ab2**2*aambb*y*dddi
1028ddxc2(3,3,3)=u2*ab2*(aapbb*ddi+ab2*e*dddi)
1029end subroutine dctocdd
1030
1039subroutine plrot(rot3,n,x,y,z)! [plrot]
1040implicit none
1041integer, intent(IN ):: n
1042real(sp),dimension(3,3),intent(IN ):: rot3
1043real(sp),dimension(n), intent(INOUT):: x,y,z
1044real(sp),dimension(3) :: t
1045integer :: k
1046do k=1,n
1047 t(1)=x(k); t(2)=y(k); t(3)=z(k)
1048 t=matmul(rot3,t)
1049 x(k)=t(1); y(k)=t(2); z(k)=t(3)
1050enddo
1051end subroutine plrot
1052
1061subroutine plroti(rot3,n,x,y,z)! [plroti]
1062implicit none
1063integer, intent(IN ):: n
1064real(sp),dimension(3,3),intent(IN ):: rot3
1065real(sp),dimension(n), intent(INOUT):: x,y,z
1066real(sp),dimension(3) :: t
1067integer :: k
1068do k=1,n
1069 t(1)=x(k); t(2)=y(k); t(3)=z(k)
1070 t=matmul(t,rot3)
1071 x(k)=t(1); y(k)=t(2); z(k)=t(3)
1072enddo
1073end subroutine plroti
1074
1083subroutine dplrot(rot3,n,x,y,z)! [plrot]
1084implicit none
1085integer, intent(IN ):: n
1086real(dP),dimension(3,3),intent(IN ):: rot3
1087real(dP),dimension(n), intent(INOUT):: x,y,z
1088real(dP),dimension(3) :: t
1089integer :: k
1090do k=1,n
1091 t(1)=x(k); t(2)=y(k); t(3)=z(k)
1092 t=matmul(rot3,t)
1093 x(k)=t(1); y(k)=t(2); z(k)=t(3)
1094enddo
1095end subroutine dplrot
1096
1105subroutine dplroti(rot3,n,x,y,z)! [plroti]
1106implicit none
1107integer, intent(IN ):: n
1108real(dP),dimension(3,3),intent(IN ):: rot3
1109real(dP),dimension(n), intent(INOUT):: x,y,z
1110real(dP),dimension(3) :: t
1111integer :: k
1112do k=1,n
1113 t(1)=x(k); t(2)=y(k); t(3)=z(k)
1114 t=matmul(t,rot3)
1115 x(k)=t(1); y(k)=t(2); z(k)=t(3)
1116enddo
1117end subroutine dplroti
1118
1127subroutine plctoc(s,n,x,y,z)! [plctoc]
1128use pietc_s, only: u1
1129implicit none
1130integer, intent(IN ):: n
1131real(sp), intent(IN ):: s
1132real(sp),dimension(n),intent(INOUT):: x,y,z
1133real(sp) :: a,b,d,e,ab2,aa,bb,di,aapbb,aambb
1134integer :: i
1135a=s+u1
1136b=s-u1
1137ab2=a*b*2
1138aa=a*a
1139bb=b*b
1140aapbb=aa+bb
1141aambb=aa-bb
1142do i=1,n
1143 d=aapbb-ab2*z(i)
1144 e=aapbb*z(i)-ab2
1145 di=1/d
1146 x(i)=(aambb*x(i))*di
1147 y(i)=(aambb*y(i))*di
1148 z(i)=e*di
1149enddo
1150end subroutine plctoc
1151
1152end module pmat5
1153
1154
1155
Constants for orientation and stretching of map.
Definition pmat5.f90:7
real(sp) sc
Schmidt stretching factor.
Definition pmat5.f90:11
real(sp), dimension(3, 3) rotm
Orthogonal rotation matrix.
Definition pmat5.f90:10
real(sp) sci
Schmidt inverse stretching factor.
Definition pmat5.f90:12
Constants for orientation and stretching of map.
Definition pmat5.f90:17
real(dp) sci
Schmidt inverse stretching factor.
Definition pmat5.f90:22
real(dp) sc
Schmidt stretching factor.
Definition pmat5.f90:21
real(dp), dimension(3, 3) rotm
Orthogonal rotation matrix.
Definition pmat5.f90:20
Some of the commonly used constants (pi etc) mainly for double-precision subroutines.
Definition pietc.f90:14
real(dp), parameter dtor
Degrees to radians.
Definition pietc.f90:54
real(dp), parameter u1
one
Definition pietc.f90:20
real(dp), parameter rtod
radians to degrees
Definition pietc.f90:55
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 sp
Single precision real kind.
Definition pkind.f90:10
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
Utility routines for orienting the globe and basic geographical mappings.
Definition pmat5.f90:27
subroutine sgtocd(dlat, dlon, xe, dxedlat, dxedlon)
Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to lati...
Definition pmat5.f90:470
subroutine sinivmap(alon0, alat0, rot3)
Initialize the rotation matrix ROT3 needed to transform standard earth-centered cartesian components ...
Definition pmat5.f90:108
subroutine dininmap(alon0, alat0, rot3)
Initialize the rotation matrix ROT3 needed to transform standard earth-centered cartesian components ...
Definition pmat5.f90:86
subroutine dctocd(s, xc1, xc2, dxc2)
Evaluate Schmidt transformation, xc1 --> xc2, with scaling parameter s, and its jacobian,...
Definition pmat5.f90:947
subroutine dgrtocd(rlat, rlon, xe, dxedlat, dxedlon)
Transform "Geographical" to "Cartesian" coordinates, together with the partial derivatives of cartesi...
Definition pmat5.f90:280
subroutine dctocdd(s, xc1, xc2, dxc2, ddxc2)
Evaluate Schmidt transformation, xc1 --> xc2, with scaling parameter s, its jacobian,...
Definition pmat5.f90:988
subroutine gtoframev(plat, plon, xp, yp, zp)
Given a geographical point by its degrees lat and lon, plat and plon, return its standard orthogonal ...
Definition pmat5.f90:662
subroutine dplroti(rot3, n, x, y, z)
Invert the rotation of a three-dimensional polyline.
Definition pmat5.f90:1106
subroutine sframetwist(sxp, syp, szp, sxv, syv, szv, stwist)
Given a principal cartesian orthonormal frame, {xp,yp,zp} (i.e., at P with Earth-centered cartesians,...
Definition pmat5.f90:751
subroutine dgtocdd(dlat, dlon, xe, dxedlat, dxedlon, ddxedlatdlat, ddxedlatdlon, ddxedlondlon)
Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to lati...
Definition pmat5.f90:565
subroutine sgrtoc(rlat, rlon, xe)
Transform "Geographical" to "Cartesian" coordinates.
Definition pmat5.f90:220
subroutine dgtoc(dlat, dlon, xe)
Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to lati...
Definition pmat5.f90:445
subroutine dgtocd(dlat, dlon, xe, dxedlat, dxedlon)
Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to lati...
Definition pmat5.f90:496
subroutine sgrtocdd(rlat, rlon, xe, dxedlat, dxedlon, ddxedlatdlat, ddxedlatdlon, ddxedlondlon)
Transform "Geographical" to "Cartesian" coordinates, together with the 1st and 2nd partial derivative...
Definition pmat5.f90:307
subroutine dctog(xe, dlat, dlon)
Transform "Cartesian" to "Geographical" coordinates, where the geographical coordinates refer to lati...
Definition pmat5.f90:398
subroutine sctog(xe, dlat, dlon)
Transform "Cartesian" to "Geographical" coordinates, where the geographical coordinates refer to lati...
Definition pmat5.f90:373
subroutine sgtoc(dlat, dlon, xe)
Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to lati...
Definition pmat5.f90:423
subroutine sparaframe(sxp, syp, szp, sxv, syv, szv)
Take a principal reference orthonormal frame, {xp,yp,zp} and a dependent point defined by unit vector...
Definition pmat5.f90:696
subroutine sctocd(s, xc1, xc2, dxc2)
Evaluate Schmidt transformation, xc1 --> xc2, with scaling parameter s, and its jacobian,...
Definition pmat5.f90:825
subroutine sctogr(xe, rlat, rlon)
Transform "Cartesian" to "Geographical" coordinates, where the geographical coordinates refer to lati...
Definition pmat5.f90:173
subroutine dctoc(s, xc1, xc2)
Evaluate Schmidt transformation, xc1 --> xc2, with scaling parameter s.
Definition pmat5.f90:916
subroutine sctoc(s, xc1, xc2)
Evaluate Schmidt transformation, xc1 --> xc2, with scaling parameter s.
Definition pmat5.f90:795
subroutine dgrtocdd(rlat, rlon, xe, dxedlat, dxedlon, ddxedlatdlat, ddxedlatdlon, ddxedlondlon)
Transform "Geographical" to "Cartesian" coordinates, together with the 1st and 2nd partial derivative...
Definition pmat5.f90:340
subroutine dgrtoc(rlat, rlon, xe)
Transform "Geographical" to "Cartesian" coordinates.
Definition pmat5.f90:237
subroutine gtoframem(plat, plon, orth)
From the degree lat and lon (plat and plon) return the standard orthogonal 3D frame at this location ...
Definition pmat5.f90:617
subroutine sininmap(alon0, alat0, rot3)
Initialize the rotation matrix ROT3 needed to transform standard earth-centered cartesian components ...
Definition pmat5.f90:64
subroutine sgtoframev(splat, splon, sxp, syp, szp)
Given a geographical point by its degrees lat and lon, plat and plon, return its standard orthogonal ...
Definition pmat5.f90:637
subroutine sctocdd(s, xc1, xc2, dxc2, ddxc2)
Evaluate Schmidt transformation, xc1 --> xc2, with scaling parameter s, its jacobian,...
Definition pmat5.f90:866
subroutine dplrot(rot3, n, x, y, z)
Apply a constant rotation to a three dimensional polyline.
Definition pmat5.f90:1084
subroutine dinivmap(alon0, alat0, rot3)
Initialize the rotation matrix ROT3 needed to transform standard earth-centered cartesian components ...
Definition pmat5.f90:140
subroutine dctogr(xe, rlat, rlon)
Transform "Cartesian" to "Geographical" coordinates, where the geographical coordinates refer to lati...
Definition pmat5.f90:198
subroutine sgrtocd(rlat, rlon, xe, dxedlat, dxedlon)
Transform "Geographical" to "Cartesian" coordinates, together with the partial derivatives of cartesi...
Definition pmat5.f90:257
subroutine sgtoframem(splat, splon, sorth)
From the degree lat and lon (plat and plon) return the standard orthogonal 3D frame at this location ...
Definition pmat5.f90:600
subroutine sgtocdd(dlat, dlon, xe, dxedlat, dxedlon, ddxedlatdlat, ddxedlatdlon, ddxedlondlon)
Transform "Geographical" to "Cartesian" coordinates, where the geographical coordinates refer to lati...
Definition pmat5.f90:528