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]

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 137 of file scale_mapproj.F90.

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

Referenced by scale_grid_real::real_setup().

137  use scale_process, only: &
139  use scale_const, only: &
140  pi_rp => const_pi, &
141  d2r_rp => const_d2r, &
142  radius_rp => const_radius
143  implicit none
144 
145  real(RP), intent(in) :: domain_center_x
146  real(RP), intent(in) :: domain_center_y
147 
148  namelist / param_mapproj / &
149  mprj_basepoint_lon, &
150  mprj_basepoint_lat, &
151  mprj_basepoint_x, &
152  mprj_basepoint_y, &
153  mprj_type, &
154  mprj_rotation, &
155  mprj_lc_lat1, &
156  mprj_lc_lat2, &
157  mprj_ps_lat, &
158  mprj_m_lat, &
159  mprj_ec_lat
160 
161  integer :: ierr
162  !---------------------------------------------------------------------------
163 
164  if( io_l ) write(io_fid_log,*)
165  if( io_l ) write(io_fid_log,*) '+++ Module[MAPPROJ]/Categ[GRID]'
166 
167  pi = real(pi_rp, kind=dp)
168  d2r = real(d2r_rp, kind=dp)
169  radius = real(radius_rp, kind=dp)
170 
171  mprj_basepoint_x = undef
172  mprj_basepoint_y = undef
173  mprj_ps_lat = undef
174  mprj_m_lat = undef
175  mprj_ec_lat = undef
176 
177  mprj_basepoint_x = domain_center_x
178  mprj_basepoint_y = domain_center_y
179 
180  !--- read namelist
181  rewind(io_fid_conf)
182  read(io_fid_conf,nml=param_mapproj,iostat=ierr)
183 
184  if( ierr < 0 ) then !--- missing
185  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
186  elseif( ierr > 0 ) then !--- fatal error
187  write(*,*) 'xxx Not appropriate names in namelist PARAM_MAPPROJ. Check!'
188  call prc_mpistop
189  endif
190  if( io_l ) write(io_fid_log,nml=param_mapproj)
191 
192 
193  if( mprj_ps_lat == undef ) mprj_ps_lat = mprj_basepoint_lat
194  if( mprj_m_lat == undef ) mprj_m_lat = mprj_basepoint_lat
195  if( mprj_ec_lat == undef ) mprj_ec_lat = mprj_basepoint_lat
196 
197  if( io_l ) write(io_fid_log,*)
198  if( io_l ) write(io_fid_log,*) '+++ MPRJ_basepoint_lon:', mprj_basepoint_lon
199  if( io_l ) write(io_fid_log,*) '+++ MPRJ_basepoint_lat:', mprj_basepoint_lat
200  if( io_l ) write(io_fid_log,*) '+++ MPRJ_basepoint_x :', mprj_basepoint_x
201  if( io_l ) write(io_fid_log,*) '+++ MPRJ_basepoint_y :', mprj_basepoint_y
202 
203  if( io_l ) write(io_fid_log,*)
204  if( io_l ) write(io_fid_log,*) '*** Projection type:', trim(mprj_type)
205  select case(trim(mprj_type))
206  case('NONE')
207  if( io_l ) write(io_fid_log,*) '=> NO map projection'
208  call mprj_none_setup
209  case('LC')
210  if( io_l ) write(io_fid_log,*) '=> Lambert Conformal projection'
211  call mprj_lambertconformal_setup
212  case('PS')
213  if( io_l ) write(io_fid_log,*) '=> Polar Stereographic projection'
214  call mprj_polarstereographic_setup
215  case('MER')
216  if( io_l ) write(io_fid_log,*) '=> Mercator projection'
217  call mprj_mercator_setup
218  case('EC')
219  if( io_l ) write(io_fid_log,*) '=> Equidistant Cylindrical projection'
220  call mprj_equidistantcylindrical_setup
221  case default
222  write(*,*) ' xxx Unsupported TYPE. STOP'
223  call prc_mpistop
224  endselect
225 
226  return
subroutine, public prc_mpistop
Abort MPI.
real(rp), public const_radius
radius of the planet [m]
Definition: scale_const.F90:46
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:35
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
real(rp), public const_pi
pi
Definition: scale_const.F90:34
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 236 of file scale_mapproj.F90.

References scale_process::prc_mpistop().

Referenced by scale_grid_real::real_update_z().

236  use scale_process, only: &
238  implicit none
239 
240  real(RP), intent(in) :: x
241  real(RP), intent(in) :: y
242  real(RP), intent(out) :: lon ! [rad]
243  real(RP), intent(out) :: lat ! [rad]
244  !---------------------------------------------------------------------------
245 
246  select case(mprj_type)
247  case('NONE')
248  call mprj_none_xy2lonlat( x, y, lon, lat )
249  case('LC')
250  call mprj_lambertconformal_xy2lonlat( x, y, lon, lat )
251  case('PS')
252  call mprj_polarstereographic_xy2lonlat( x, y, lon, lat )
253  case('MER')
254  call mprj_mercator_xy2lonlat( x, y, lon, lat )
255  case('EC')
256  call mprj_equidistantcylindrical_xy2lonlat( x, y, lon, lat )
257  case default
258  write(*,*) ' xxx Unsupported TYPE. STOP'
259  call prc_mpistop
260  endselect
261 
262  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 272 of file scale_mapproj.F90.

References scale_process::prc_mpistop().

272  use scale_process, only: &
274  implicit none
275 
276  real(RP), intent(in) :: lon ! [rad]
277  real(RP), intent(in) :: lat ! [rad]
278  real(RP), intent(out) :: x
279  real(RP), intent(out) :: y
280  !---------------------------------------------------------------------------
281 
282  select case(mprj_type)
283  case('NONE')
284  call mprj_none_lonlat2xy( lon, lat, x, y )
285  case('LC')
286  call mprj_lambertconformal_lonlat2xy( lon, lat, x, y )
287  case('PS')
288  call mprj_polarstereographic_lonlat2xy( lon, lat, x, y )
289  case('MER')
290  call mprj_mercator_lonlat2xy( lon, lat, x, y )
291  case('EC')
292  call mprj_equidistantcylindrical_lonlat2xy( lon, lat, x, y )
293  case default
294  write(*,*) ' xxx Unsupported TYPE. STOP'
295  call prc_mpistop
296  endselect
297 
298  return
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
Here is the call 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 307 of file scale_mapproj.F90.

References scale_process::prc_mpistop().

Referenced by scale_gridtrans::gtrans_setup().

307  use scale_process, only: &
309  implicit none
310 
311  real(RP), intent(in) :: lat(ia,ja) ! [rad]
312  real(RP), intent(out) :: m1 (ia,ja)
313  real(RP), intent(out) :: m2 (ia,ja)
314  !---------------------------------------------------------------------------
315 
316  select case(mprj_type)
317  case('NONE')
318  call mprj_none_mapfactor( lat, m1, m2 )
319  case('LC')
320  call mprj_lambertconformal_mapfactor( lat, m1, m2 )
321  case('PS')
322  call mprj_polarstereographic_mapfactor( lat, m1, m2 )
323  case('MER')
324  call mprj_mercator_mapfactor( lat, m1, m2 )
325  case('EC')
326  call mprj_equidistantcylindrical_mapfactor( lat, m1, m2 )
327  case default
328  write(*,*) ' xxx Unsupported TYPE. STOP'
329  call prc_mpistop
330  endselect
331 
332  return
subroutine, public prc_mpistop
Abort MPI.
integer, public ia
of x whole cells (local, with HALO)
module PROCESS
integer, public ja
of y whole cells (local, with HALO)
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 342 of file scale_mapproj.F90.

References scale_process::prc_mpistop().

342  use scale_process, only: &
344  implicit none
345 
346  real(RP), intent(out) :: rotc(2)
347  real(RP), intent(in) :: lon ! [rad]
348  real(RP), intent(in) :: lat ! [rad]
349  !---------------------------------------------------------------------------
350 
351  select case(mprj_type)
352  case('NONE')
353  call mprj_none_rotcoef_0d( rotc )
354  case('LC')
355  call mprj_lambertconformal_rotcoef_0d( rotc, lon, lat )
356  case('PS')
357  call mprj_polarstereographic_rotcoef_0d( rotc, lon, lat )
358  case('MER')
359  call mprj_mercator_rotcoef_0d( rotc )
360  case('EC')
361  call mprj_equidistantcylindrical_rotcoef_0d( rotc )
362  case default
363  write(*,*) ' xxx Unsupported TYPE. STOP'
364  call prc_mpistop
365  endselect
366 
367  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 377 of file scale_mapproj.F90.

References scale_grid_index::ia, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::ja, dc_log::log(), mprj_basepoint_lat, mprj_basepoint_lon, and scale_process::prc_mpistop().

377  use scale_process, only: &
379  implicit none
380 
381  real(RP), intent(out) :: rotc(ia,ja,2)
382  real(RP), intent(in) :: lon (ia,ja) ! [rad]
383  real(RP), intent(in) :: lat (ia,ja) ! [rad]
384  !---------------------------------------------------------------------------
385 
386  select case(mprj_type)
387  case('NONE')
388  call mprj_none_rotcoef_2d( rotc )
389  case('LC')
390  call mprj_lambertconformal_rotcoef_2d( rotc, lon, lat )
391  case('PS')
392  call mprj_polarstereographic_rotcoef_2d( rotc, lon, lat )
393  case('MER')
394  call mprj_mercator_rotcoef_2d( rotc )
395  case('EC')
396  call mprj_equidistantcylindrical_rotcoef_2d( rotc )
397  case default
398  write(*,*) ' xxx Unsupported TYPE. STOP'
399  call prc_mpistop
400  endselect
401 
402  return
subroutine, public prc_mpistop
Abort MPI.
integer, public ia
of x whole cells (local, with HALO)
module PROCESS
integer, public ja
of y whole cells (local, with HALO)
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 48 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().

48  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 49 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().

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