SCALE-RM
Data Types | Functions/Subroutines | Variables
scale_mapproj Module Reference

module Map projection More...

Functions/Subroutines

subroutine, public mprj_setup (DOMAIN_CENTER_X, DOMAIN_CENTER_Y)
 Setup. More...
 
subroutine, public mprj_xy2lonlat (x, y, lon, lat)
 (x,y) -> (lon,lat) More...
 
subroutine, public mprj_lonlat2xy (lon, lat, x, y)
 (lon,lat) -> (x,y) More...
 
subroutine, public mprj_mapfactor (lat, m1, m2)
 (x,y) -> (lon,lat) More...
 
subroutine mprj_rotcoef_0d (rotc, lon, lat)
 u(lat,lon) = cos u(x,y) - sin v(x,y) v(lat,lon) = sin u(x,y) + cos v(x,y) More...
 
subroutine mprj_rotcoef_2d (rotc, lon, lat)
 u(lat,lon) = cos u(x,y) - sin v(x,y) v(lat,lon) = sin u(x,y) + cos v(x,y) More...
 

Variables

real(rp), public mprj_basepoint_lon = 135.221_RP
 
real(rp), public mprj_basepoint_lat = 34.653_RP
 

Detailed Description

module Map projection

Description
Map projection module
Author
Team SCALE
History
  • 2013-10-24 (H.Yashiro) [new]
NAMELIST
  • PARAM_MAPPROJ
    nametypedefault valuecomment
    MPRJ_BASEPOINT_LON real(RP) 135.221_RP position of base point (domain center) in real world [deg]
    MPRJ_BASEPOINT_LAT real(RP) 34.653_RP position of base point (domain center) in real world [deg]
    MPRJ_BASEPOINT_X real(DP) position of base point in the model [m]
    MPRJ_BASEPOINT_Y real(DP) position of base point in the model [m]
    MPRJ_TYPE character(len=H_SHORT) 'NONE' map projection type
    MPRJ_ROTATION real(DP) 0.0_DP rotation factor (only for 'NONE' type)
    MPRJ_LC_LAT1 real(DP) 30.0_DP standard latitude1 for L.C. projection [deg]
    MPRJ_LC_LAT2 real(DP) 60.0_DP standard latitude2 for L.C. projection [deg]
    MPRJ_PS_LAT real(DP) standard latitude1 for P.S. projection [deg]
    MPRJ_M_LAT real(DP) 0.0_DP standard latitude1 for Mer. projection [deg]
    MPRJ_EC_LAT real(DP) 0.0_DP standard latitude1 for E.C. projection [deg]

History Output
No history output

Function/Subroutine Documentation

◆ mprj_setup()

subroutine, public scale_mapproj::mprj_setup ( real(rp), intent(in)  DOMAIN_CENTER_X,
real(rp), intent(in)  DOMAIN_CENTER_Y 
)

Setup.

Parameters
[in]domain_center_xcenter position of global domain [m]: x
[in]domain_center_ycenter position of global domain [m]: y

Definition at line 134 of file scale_mapproj.F90.

References scale_const::const_d2r, scale_const::const_pi, scale_const::const_radius, scale_const::const_undef, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_fid_nml, scale_stdio::io_l, scale_stdio::io_nml, mprj_basepoint_lat, mprj_basepoint_lon, and scale_process::prc_mpistop().

Referenced by scale_grid_real::real_setup().

134  use scale_process, only: &
136  use scale_const, only: &
137  undef => const_undef, &
138  pi_rp => const_pi, &
139  d2r_rp => const_d2r, &
140  radius_rp => const_radius
141  implicit none
142 
143  real(RP), intent(in) :: DOMAIN_CENTER_X
144  real(RP), intent(in) :: DOMAIN_CENTER_Y
145 
146  namelist / param_mapproj / &
147  mprj_basepoint_lon, &
148  mprj_basepoint_lat, &
149  mprj_basepoint_x, &
150  mprj_basepoint_y, &
151  mprj_type, &
152  mprj_rotation, &
153  mprj_lc_lat1, &
154  mprj_lc_lat2, &
155  mprj_ps_lat, &
156  mprj_m_lat, &
157  mprj_ec_lat
158 
159  integer :: ierr
160  !---------------------------------------------------------------------------
161 
162  if( io_l ) write(io_fid_log,*)
163  if( io_l ) write(io_fid_log,*) '++++++ Module[MAPPROJ] / Categ[ATMOS-RM GRID] / Origin[SCALElib]'
164 
165  pi = real(PI_RP, kind=dp)
166  d2r = real(D2R_RP, kind=dp)
167  radius = real(RADIUS_RP, kind=dp)
168 
169  mprj_basepoint_x = undef
170  mprj_basepoint_y = undef
171  mprj_ps_lat = undef
172 
173  mprj_basepoint_x = domain_center_x
174  mprj_basepoint_y = domain_center_y
175 
176  !--- read namelist
177  rewind(io_fid_conf)
178  read(io_fid_conf,nml=param_mapproj,iostat=ierr)
179 
180  if( ierr < 0 ) then !--- missing
181  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
182  elseif( ierr > 0 ) then !--- fatal error
183  write(*,*) 'xxx Not appropriate names in namelist PARAM_MAPPROJ. Check!'
184  call prc_mpistop
185  endif
186  if( io_nml ) write(io_fid_nml,nml=param_mapproj)
187 
188  if( mprj_ps_lat == undef ) mprj_ps_lat = mprj_basepoint_lat
189 
190  if( io_l ) write(io_fid_log,*)
191  if( io_l ) write(io_fid_log,*) '*** Map projection type : ', trim(mprj_type)
192  select case(trim(mprj_type))
193  case('NONE')
194  if( io_l ) write(io_fid_log,*) '*** => NO map projection'
195  call mprj_none_setup
196  case('LC')
197  if( io_l ) write(io_fid_log,*) '*** => Lambert Conformal projection'
198  call mprj_lambertconformal_setup
199  case('PS')
200  if( io_l ) write(io_fid_log,*) '*** => Polar Stereographic projection'
201  call mprj_polarstereographic_setup
202  case('MER')
203  if( io_l ) write(io_fid_log,*) '*** => Mercator projection'
204  call mprj_mercator_setup
205  case('EC')
206  if( io_l ) write(io_fid_log,*) '*** => Equidistant Cylindrical projection'
207  call mprj_equidistantcylindrical_setup
208  case default
209  write(*,*) 'xxx Unsupported MPRJ_type. STOP'
210  call prc_mpistop
211  endselect
212 
213  if( io_l ) write(io_fid_log,'(1x,A,F15.3)') '*** Basepoint(x) = ', mprj_basepoint_x
214  if( io_l ) write(io_fid_log,'(1x,A,F15.3)') '*** Basepoint(y) = ', mprj_basepoint_y
215 
216  return
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, parameter, public dp
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mprj_xy2lonlat()

subroutine, public scale_mapproj::mprj_xy2lonlat ( real(rp), intent(in)  x,
real(rp), intent(in)  y,
real(rp), intent(out)  lon,
real(rp), intent(out)  lat 
)

(x,y) -> (lon,lat)

Definition at line 226 of file scale_mapproj.F90.

References scale_process::prc_mpistop().

Referenced by scale_grid_real::real_update_z().

226  use scale_process, only: &
228  implicit none
229 
230  real(RP), intent(in) :: x
231  real(RP), intent(in) :: y
232  real(RP), intent(out) :: lon ! [rad]
233  real(RP), intent(out) :: lat ! [rad]
234  !---------------------------------------------------------------------------
235 
236  select case(mprj_type)
237  case('NONE')
238  call mprj_none_xy2lonlat( x, y, lon, lat )
239  case('LC')
240  call mprj_lambertconformal_xy2lonlat( x, y, lon, lat )
241  case('PS')
242  call mprj_polarstereographic_xy2lonlat( x, y, lon, lat )
243  case('MER')
244  call mprj_mercator_xy2lonlat( x, y, lon, lat )
245  case('EC')
246  call mprj_equidistantcylindrical_xy2lonlat( x, y, lon, lat )
247  case default
248  write(*,*) 'xxx Unsupported MPRJ_type. STOP'
249  call prc_mpistop
250  endselect
251 
252  return
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mprj_lonlat2xy()

subroutine, public scale_mapproj::mprj_lonlat2xy ( real(rp), intent(in)  lon,
real(rp), intent(in)  lat,
real(rp), intent(out)  x,
real(rp), intent(out)  y 
)

(lon,lat) -> (x,y)

Definition at line 262 of file scale_mapproj.F90.

References scale_process::prc_mpistop().

Referenced by scale_gridtrans::gtrans_rotcoef().

262  use scale_process, only: &
264  implicit none
265 
266  real(RP), intent(in) :: lon ! [rad]
267  real(RP), intent(in) :: lat ! [rad]
268  real(RP), intent(out) :: x
269  real(RP), intent(out) :: y
270  !---------------------------------------------------------------------------
271 
272  select case(mprj_type)
273  case('NONE')
274  call mprj_none_lonlat2xy( lon, lat, x, y )
275  case('LC')
276  call mprj_lambertconformal_lonlat2xy( lon, lat, x, y )
277  case('PS')
278  call mprj_polarstereographic_lonlat2xy( lon, lat, x, y )
279  case('MER')
280  call mprj_mercator_lonlat2xy( lon, lat, x, y )
281  case('EC')
282  call mprj_equidistantcylindrical_lonlat2xy( lon, lat, x, y )
283  case default
284  write(*,*) 'xxx Unsupported MPRJ_type. STOP'
285  call prc_mpistop
286  endselect
287 
288  return
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mprj_mapfactor()

subroutine, public scale_mapproj::mprj_mapfactor ( real(rp), dimension(ia,ja), intent(in)  lat,
real(rp), dimension (ia,ja), intent(out)  m1,
real(rp), dimension (ia,ja), intent(out)  m2 
)

(x,y) -> (lon,lat)

Definition at line 297 of file scale_mapproj.F90.

References scale_process::prc_mpistop().

Referenced by scale_gridtrans::gtrans_setup().

297  use scale_process, only: &
299  implicit none
300 
301  real(RP), intent(in) :: lat(IA,JA) ! [rad]
302  real(RP), intent(out) :: m1 (IA,JA)
303  real(RP), intent(out) :: m2 (IA,JA)
304  !---------------------------------------------------------------------------
305 
306  select case(mprj_type)
307  case('NONE')
308  call mprj_none_mapfactor( lat, m1, m2 )
309  case('LC')
310  call mprj_lambertconformal_mapfactor( lat, m1, m2 )
311  case('PS')
312  call mprj_polarstereographic_mapfactor( lat, m1, m2 )
313  case('MER')
314  call mprj_mercator_mapfactor( lat, m1, m2 )
315  case('EC')
316  call mprj_equidistantcylindrical_mapfactor( lat, m1, m2 )
317  case default
318  write(*,*) 'xxx Unsupported MPRJ_type. STOP'
319  call prc_mpistop
320  endselect
321 
322  return
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mprj_rotcoef_0d()

subroutine scale_mapproj::mprj_rotcoef_0d ( real(rp), dimension(2), intent(out)  rotc,
real(rp), intent(in)  lon,
real(rp), intent(in)  lat 
)

u(lat,lon) = cos u(x,y) - sin v(x,y) v(lat,lon) = sin u(x,y) + cos v(x,y)

Parameters
[out]rotcrotc(:,:,1)->cos, rotc(:,:,2)->sin

Definition at line 332 of file scale_mapproj.F90.

References scale_process::prc_mpistop().

332  use scale_process, only: &
334  implicit none
335 
336  real(RP), intent(out) :: rotc(2)
337  real(RP), intent(in) :: lon ! [rad]
338  real(RP), intent(in) :: lat ! [rad]
339  !---------------------------------------------------------------------------
340 
341  select case(mprj_type)
342  case('NONE')
343  call mprj_none_rotcoef_0d( rotc )
344  case('LC')
345  call mprj_lambertconformal_rotcoef_0d( rotc, lon, lat )
346  case('PS')
347  call mprj_polarstereographic_rotcoef_0d( rotc, lon, lat )
348  case('MER')
349  call mprj_mercator_rotcoef_0d( rotc )
350  case('EC')
351  call mprj_equidistantcylindrical_rotcoef_0d( rotc )
352  case default
353  write(*,*) 'xxx Unsupported MPRJ_type. STOP'
354  call prc_mpistop
355  endselect
356 
357  return
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
Here is the call graph for this function:

◆ mprj_rotcoef_2d()

subroutine scale_mapproj::mprj_rotcoef_2d ( real(rp), dimension(ia,ja,2), intent(out)  rotc,
real(rp), dimension (ia,ja), intent(in)  lon,
real(rp), dimension (ia,ja), intent(in)  lat 
)

u(lat,lon) = cos u(x,y) - sin v(x,y) v(lat,lon) = sin u(x,y) + cos v(x,y)

Parameters
[out]rotcrotc(:,:,1)->cos, rotc(:,:,2)->sin

Definition at line 367 of file scale_mapproj.F90.

References scale_stdio::io_fid_log, scale_stdio::io_l, mprj_basepoint_lat, mprj_basepoint_lon, and scale_process::prc_mpistop().

367  use scale_process, only: &
369  implicit none
370 
371  real(RP), intent(out) :: rotc(IA,JA,2)
372  real(RP), intent(in) :: lon (IA,JA) ! [rad]
373  real(RP), intent(in) :: lat (IA,JA) ! [rad]
374  !---------------------------------------------------------------------------
375 
376  select case(mprj_type)
377  case('NONE')
378  call mprj_none_rotcoef_2d( rotc )
379  case('LC')
380  call mprj_lambertconformal_rotcoef_2d( rotc, lon, lat )
381  case('PS')
382  call mprj_polarstereographic_rotcoef_2d( rotc, lon, lat )
383  case('MER')
384  call mprj_mercator_rotcoef_2d( rotc )
385  case('EC')
386  call mprj_equidistantcylindrical_rotcoef_2d( rotc )
387  case default
388  write(*,*) 'xxx Unsupported MPRJ_type. STOP'
389  call prc_mpistop
390  endselect
391 
392  return
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
Here is the call graph for this function:

Variable Documentation

◆ mprj_basepoint_lon

real(rp), public scale_mapproj::mprj_basepoint_lon = 135.221_RP

Definition at line 45 of file scale_mapproj.F90.

Referenced by mprj_rotcoef_2d(), mprj_setup(), scale_grid_real::real_update_z(), and mod_urban_phy_driver::urban_phy_driver().

45  real(RP), public :: MPRJ_basepoint_lon = 135.221_rp ! position of base point (domain center) in real world [deg]

◆ mprj_basepoint_lat

real(rp), public scale_mapproj::mprj_basepoint_lat = 34.653_RP

Definition at line 46 of file scale_mapproj.F90.

Referenced by mprj_rotcoef_2d(), mprj_setup(), scale_grid_real::real_update_z(), and mod_urban_phy_driver::urban_phy_driver().

46  real(RP), public :: MPRJ_basepoint_lat = 34.653_rp ! position of base point (domain center) in real world [deg]