SICOPOLIS V3.1
 All Classes Files Functions Variables Macros
num_coord.F90
Go to the documentation of this file.
1 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 !
3 ! Subroutine : n u m _ c o o r d
4 !
5 !> @file
6 !!
7 !! Computation of position (x,y) in the numerical domain for longitude lambda
8 !! and latitude phi (in rad).
9 !!
10 !! @section Copyright
11 !!
12 !! Copyright 2009-2013 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 position (x,y) in the numerical domain for longitude lambda
35 !! and latitude phi (in rad).
36 !<------------------------------------------------------------------------------
37 subroutine num_coord(lambda_val, phi_val, x_val, y_val)
38 
39 use sico_types
41 
42 implicit none
43 
44 real(dp), intent(in) :: lambda_val, phi_val
45 
46 real(dp), intent(out) :: x_val, y_val
47 
48 #if (GRID==0 || GRID==1)
49 
50 !-------- Stereographic projection --------
51 
52 call stereo_forw_ellipsoid(lambda_val, phi_val, a, b, &
53  lambda0, phi0, x_val, y_val)
54 
55 #elif GRID==2
56 
57 !-------- Identity --------
58 
59 x_val = lambda_val
60 y_val = phi_val
61 
62 #endif
63 
64 end subroutine num_coord
65 
66 !-------------------------------------------------------------------------------
67 !> Forward stereographic projection for an ellipsoidal planet.
68 !<------------------------------------------------------------------------------
69 subroutine stereo_forw_ellipsoid(lambda_val, phi_val, A, B, &
70  lambda0, phi0, x_val, y_val)
71 
72 use sico_types
73 
74 implicit none
75 real(dp), intent(in) :: lambda_val, phi_val, a, b, lambda0, phi0
76 real(dp), intent(out) :: x_val, y_val
77 
78 integer(i4b) :: l
79 integer(i4b) :: sign_phi0
80 real(dp) :: phi_aux, phi0_aux
81 real(dp) :: e, mc, t, tc, kp, rho, phi_p
82 real(dp) :: sinphi0, sinlambda0, cosphi0, coslambda0
83 
84 real(dp), parameter :: pi = 3.141592653589793_dp, &
85  eps = 1.0e-05_dp
86 
87 if (phi0 > eps) then ! for northern hemisphere
88  sign_phi0 = 1
89 else if (phi0 < (-eps)) then ! for southern hemisphere
90  sign_phi0 = -1
91 else
92  stop ' num_coord: phi0 must be different from zero!'
93 end if
94 
95 phi_aux = phi_val * sign_phi0
96 phi0_aux = phi0 * sign_phi0
97 
98 e=sqrt((a**2-b**2)/(a**2))
99 
100 sinphi0 = sin(phi0_aux)
101 sinlambda0 = sin(lambda0)
102 cosphi0 = cos(phi0_aux)
103 coslambda0 = cos(lambda0)
104 
105 mc=cosphi0/sqrt(1.0_dp-e*e*sinphi0*sinphi0)
106 t=sqrt(((1.0_dp-sin(phi_aux))/(1.0_dp+sin(phi_aux)))* &
107  ((1.0_dp+e*sin(phi_aux))/(1.0_dp-e*sin(phi_aux)))**e)
108 tc=sqrt(((1.0_dp-sinphi0)/(1.0_dp+sinphi0))* &
109  ((1.0_dp+e*sinphi0)/(1.0_dp-e*sinphi0))**e)
110 rho=a*mc*t/tc
111 
112 x_val = rho*sin(lambda_val-lambda0)
113 y_val = -sign_phi0 * rho*cos(lambda_val-lambda0)
114 
115 end subroutine stereo_forw_ellipsoid
116 !