SICOPOLIS V3.2
 All Classes Files Functions Variables Macros
geo_coord.F90
Go to the documentation of this file.
1 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 !
3 ! Subroutine : g e o _ c o o r d
4 !
5 !> @file
6 !!
7 !! Computation of longitude lambda and latitude phi for position (x,y) in the
8 !! numerical domain.
9 !!
10 !! @section Copyright
11 !!
12 !! Copyright 2009-2016 Ralf Greve, Reinhard Calov, Alex Robinson
13 !!
14 !! @section License
15 !!
16 !! This file is part of SICOPOLIS.
17 !!
18 !! SICOPOLIS is free software: you can redistribute it and/or modify
19 !! it under the terms of the GNU General Public License as published by
20 !! the Free Software Foundation, either version 3 of the License, or
21 !! (at your option) any later version.
22 !!
23 !! SICOPOLIS is distributed in the hope that it will be useful,
24 !! but WITHOUT ANY WARRANTY; without even the implied warranty of
25 !! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
26 !! GNU General Public License for more details.
27 !!
28 !! You should have received a copy of the GNU General Public License
29 !! along with SICOPOLIS. If not, see <http://www.gnu.org/licenses/>.
30 !<
31 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
32 
33 !-------------------------------------------------------------------------------
34 !> Computation of longitude lambda and latitude phi for position (x,y) in the
35 !! numerical domain.
36 !<------------------------------------------------------------------------------
37 subroutine geo_coord(phi_val, lambda_val, x_val, y_val)
38 
39 use sico_types
41 use sico_vars
42 
43 implicit none
44 
45 real(dp), intent(in) :: x_val, y_val
46 
47 real(dp), intent(out) :: lambda_val, phi_val
48 
49 real(dp) :: k
50 
51 #if (GRID==0 || GRID==1)
52 
53 !-------- Inverse stereographic projection --------
54 
55 call stereo_inv_ellipsoid(x_val, y_val, a, b, &
56  lambda0, phi0, lambda_val, phi_val)
57 
58 #elif GRID==2
59 
60 !-------- Identity --------
61 
62 lambda_val = x_val
63 phi_val = y_val
64 
65 #endif
66 
67 end subroutine geo_coord
68 
69 !-------------------------------------------------------------------------------
70 !> Inverse stereographic projection for an ellipsoidal planet.
71 !<------------------------------------------------------------------------------
72 subroutine stereo_inv_ellipsoid(x_val, y_val, A, B, &
73  lambda0, phi0, lambda_val, phi_val)
74 
75 use sico_types
76 
77 implicit none
78 real(dp), intent(in) :: x_val, y_val, a, b, lambda0, phi0
79 real(dp), intent(out) :: lambda_val, phi_val
80 
81 integer(i4b) :: l
82 integer(i4b) :: sign_phi0
83 real(dp) :: phi_aux, phi0_aux
84 real(dp) :: e, mc, t, tc, kp, rho, phi_p, residual
85 real(dp) :: sinphi0, sinlambda0, cosphi0, coslambda0
86 
87 real(dp), parameter :: pi = 3.141592653589793_dp, &
88  eps = 1.0e-05_dp, &
89  eps_residual = 1.0e-09_dp
90 
91 if (phi0 > eps) then ! for northern hemisphere
92  sign_phi0 = 1
93 else if (phi0 < (-eps)) then ! for southern hemisphere
94  sign_phi0 = -1
95 else
96  stop ' geo_coord: phi0 must be different from zero!'
97 end if
98 
99 phi0_aux = phi0 * sign_phi0
100 
101 e=sqrt((a**2-b**2)/(a**2))
102 
103 sinphi0 = sin(phi0_aux)
104 sinlambda0 = sin(lambda0)
105 cosphi0 = cos(phi0_aux)
106 coslambda0 = cos(lambda0)
107 
108 tc=sqrt(((1.0_dp-sinphi0)/(1.0_dp+sinphi0))* &
109  ((1.0_dp+e*sinphi0)/(1.0_dp-e*sinphi0))**e)
110 mc=cosphi0/sqrt(1.0_dp-e*e*sinphi0*sinphi0)
111 rho=sqrt(x_val*x_val+y_val*y_val)
112 t=rho*tc/(a*mc)
113 
114 lambda_val = lambda0 + sign_phi0*atan2(y_val,x_val) + 0.5_dp*pi
115 
116 ! fix point iteration
117 
118 phi_p=0.5_dp*pi-2.0_dp*atan(t)
119 l=0
120 residual=3600.0_dp
121 do while(residual >= eps_residual)
122  l=l+1
123  phi_aux=0.5_dp*pi-2.0_dp*atan(t*((1.0_dp-e*sin(phi_p))/ &
124  (1.0_dp+e*sin(phi_p)))**(0.5_dp*e))
125  residual=abs(phi_aux-phi_p)
126  phi_p=phi_aux
127 end do
128 
129 phi_val = phi_aux * sign_phi0
130 
131 if (lambda_val < 0.0_dp) then
132  lambda_val = lambda_val + 2.0_dp*pi
133 else if (lambda_val >= (2.0_dp*pi)) then
134  lambda_val = lambda_val - 2.0_dp*pi
135 end if
136 
137 end subroutine stereo_inv_ellipsoid
138 !
Declarations of kind types for SICOPOLIS.
Definition: sico_types.F90:35
Declarations of global variables for SICOPOLIS (for the ANT domain).
Definition: sico_vars.F90:35
subroutine stereo_inv_ellipsoid(x_val, y_val, A, B, lambda0, phi0, lambda_val, phi_val)
Inverse stereographic projection for an ellipsoidal planet.
Definition: geo_coord.F90:72
subroutine geo_coord(phi_val, lambda_val, x_val, y_val)
Computation of longitude lambda and latitude phi for position (x,y) in the numerical domain...
Definition: geo_coord.F90:37
Declarations of global variables for SICOPOLIS.