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