PageRenderTime 39ms CodeModel.GetById 5ms RepoModel.GetById 0ms app.codeStats 0ms

/WPS/geogrid/src/module_map_utils.F

http://github.com/jbeezley/wrf-fire
FORTRAN Legacy | 2210 lines | 1165 code | 507 blank | 538 comment | 37 complexity | 9089163e9f96298a6c09e995f476fd9e MD5 | raw file
Possible License(s): AGPL-1.0
  1. MODULE map_utils
  2. ! Module that defines constants, data structures, and
  3. ! subroutines used to convert grid indices to lat/lon
  4. ! and vice versa.
  5. !
  6. ! SUPPORTED PROJECTIONS
  7. ! ---------------------
  8. ! Cylindrical Lat/Lon (code = PROJ_LATLON)
  9. ! Mercator (code = PROJ_MERC)
  10. ! Lambert Conformal (code = PROJ_LC)
  11. ! Gaussian (code = PROJ_GAUSS)
  12. ! Polar Stereographic (code = PROJ_PS)
  13. ! Rotated Lat/Lon (code = PROJ_ROTLL)
  14. !
  15. ! REMARKS
  16. ! -------
  17. ! The routines contained within were adapted from routines
  18. ! obtained from NCEP's w3 library. The original NCEP routines were less
  19. ! flexible (e.g., polar-stereo routines only supported truelat of 60N/60S)
  20. ! than what we needed, so modifications based on equations in Hoke, Hayes, and
  21. ! Renninger (AFGWC/TN/79-003) were added to improve the flexibility.
  22. ! Additionally, coding was improved to F90 standards and the routines were
  23. ! combined into this module.
  24. !
  25. ! ASSUMPTIONS
  26. ! -----------
  27. ! Grid Definition:
  28. ! For mercator, lambert conformal, and polar-stereographic projections,
  29. ! the routines within assume the following:
  30. !
  31. ! 1. Grid is dimensioned (i,j) where i is the East-West direction,
  32. ! positive toward the east, and j is the north-south direction,
  33. ! positive toward the north.
  34. ! 2. Origin is at (1,1) and is located at the southwest corner,
  35. ! regardless of hemispere.
  36. ! 3. Grid spacing (dx) is always positive.
  37. ! 4. Values of true latitudes must be positive for NH domains
  38. ! and negative for SH domains.
  39. !
  40. ! For the latlon and Gaussian projection, the grid origin may be at any
  41. ! of the corners, and the deltalat and deltalon values can be signed to
  42. ! account for this using the following convention:
  43. ! Origin Location Deltalat Sign Deltalon Sign
  44. ! --------------- ------------- -------------
  45. ! SW Corner + +
  46. ! NE Corner - -
  47. ! NW Corner - +
  48. ! SE Corner + -
  49. !
  50. ! Data Definitions:
  51. ! 1. Any arguments that are a latitude value are expressed in
  52. ! degrees north with a valid range of -90 -> 90
  53. ! 2. Any arguments that are a longitude value are expressed in
  54. ! degrees east with a valid range of -180 -> 180.
  55. ! 3. Distances are in meters and are always positive.
  56. ! 4. The standard longitude (stdlon) is defined as the longitude
  57. ! line which is parallel to the grid's y-axis (j-direction), along
  58. ! which latitude increases (NOT the absolute value of latitude, but
  59. ! the actual latitude, such that latitude increases continuously
  60. ! from the south pole to the north pole) as j increases.
  61. ! 5. One true latitude value is required for polar-stereographic and
  62. ! mercator projections, and defines at which latitude the
  63. ! grid spacing is true. For lambert conformal, two true latitude
  64. ! values must be specified, but may be set equal to each other to
  65. ! specify a tangent projection instead of a secant projection.
  66. !
  67. ! USAGE
  68. ! -----
  69. ! To use the routines in this module, the calling routines must have the
  70. ! following statement at the beginning of its declaration block:
  71. ! USE map_utils
  72. !
  73. ! The use of the module not only provides access to the necessary routines,
  74. ! but also defines a structure of TYPE (proj_info) that can be used
  75. ! to declare a variable of the same type to hold your map projection
  76. ! information. It also defines some integer parameters that contain
  77. ! the projection codes so one only has to use those variable names rather
  78. ! than remembering the acutal code when using them. The basic steps are
  79. ! as follows:
  80. !
  81. ! 1. Ensure the "USE map_utils" is in your declarations.
  82. ! 2. Declare the projection information structure as type(proj_info):
  83. ! TYPE(proj_info) :: proj
  84. ! 3. Populate your structure by calling the map_set routine:
  85. ! CALL map_set(code,lat1,lon1,knowni,knownj,dx,stdlon,truelat1,truelat2,proj)
  86. ! where:
  87. ! code (input) = one of PROJ_LATLON, PROJ_MERC, PROJ_LC, PROJ_PS,
  88. ! PROJ_GAUSS, or PROJ_ROTLL
  89. ! lat1 (input) = Latitude of grid origin point (i,j)=(1,1)
  90. ! (see assumptions!)
  91. ! lon1 (input) = Longitude of grid origin
  92. ! knowni (input) = origin point, x-location
  93. ! knownj (input) = origin point, y-location
  94. ! dx (input) = grid spacing in meters (ignored for LATLON projections)
  95. ! stdlon (input) = Standard longitude for PROJ_PS and PROJ_LC,
  96. ! deltalon (see assumptions) for PROJ_LATLON,
  97. ! ignored for PROJ_MERC
  98. ! truelat1 (input) = 1st true latitude for PROJ_PS, PROJ_LC, and
  99. ! PROJ_MERC, deltalat (see assumptions) for PROJ_LATLON
  100. ! truelat2 (input) = 2nd true latitude for PROJ_LC,
  101. ! ignored for all others.
  102. ! proj (output) = The structure of type (proj_info) that will be fully
  103. ! populated after this call
  104. !
  105. ! 4. Now that the proj structure is populated, you may call either
  106. ! of the following routines:
  107. !
  108. ! latlon_to_ij(proj, lat, lon, i, j)
  109. ! ij_to_latlon(proj, i, j, lat, lon)
  110. !
  111. ! It is incumbent upon the calling routine to determine whether or
  112. ! not the values returned are within your domain's bounds. All values
  113. ! of i, j, lat, and lon are REAL values.
  114. !
  115. !
  116. ! REFERENCES
  117. ! ----------
  118. ! Hoke, Hayes, and Renninger, "Map Preojections and Grid Systems for
  119. ! Meteorological Applications." AFGWC/TN-79/003(Rev), Air Weather
  120. ! Service, 1985.
  121. !
  122. ! NCAR MM5v3 Modeling System, REGRIDDER program, module_first_guess_map.F
  123. ! NCEP routines w3fb06, w3fb07, w3fb08, w3fb09, w3fb11, w3fb12
  124. !
  125. ! HISTORY
  126. ! -------
  127. ! 27 Mar 2001 - Original Version
  128. ! Brent L. Shaw, NOAA/FSL (CSU/CIRA)
  129. !
  130. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  131. use constants_module
  132. use misc_definitions_module
  133. use module_debug
  134. ! Define some private constants
  135. INTEGER, PRIVATE, PARAMETER :: HIGH = 8
  136. TYPE proj_info
  137. INTEGER :: code ! Integer code for projection TYPE
  138. INTEGER :: nlat ! For Gaussian -- number of latitude points
  139. ! north of the equator
  140. INTEGER :: nlon !
  141. !
  142. INTEGER :: nxmin ! Starting x-coordinate of periodic, regular lat/lon dataset
  143. INTEGER :: nxmax ! Ending x-coordinate of periodic, regular lat/lon dataset
  144. INTEGER :: ixdim ! For Rotated Lat/Lon -- number of mass points
  145. ! in an odd row
  146. INTEGER :: jydim ! For Rotated Lat/Lon -- number of rows
  147. INTEGER :: stagger ! For Rotated Lat/Lon -- mass or velocity grid
  148. REAL :: phi ! For Rotated Lat/Lon -- domain half-extent in
  149. ! degrees latitude
  150. REAL :: lambda ! For Rotated Lat/Lon -- domain half-extend in
  151. ! degrees longitude
  152. REAL :: lat1 ! SW latitude (1,1) in degrees (-90->90N)
  153. REAL :: lon1 ! SW longitude (1,1) in degrees (-180->180E)
  154. REAL :: lat0 ! For Cassini, latitude of projection pole
  155. REAL :: lon0 ! For Cassini, longitude of projection pole
  156. REAL :: dx ! Grid spacing in meters at truelats, used
  157. ! only for ps, lc, and merc projections
  158. REAL :: dy ! Grid spacing in meters at truelats, used
  159. ! only for ps, lc, and merc projections
  160. REAL :: latinc ! Latitude increment for cylindrical lat/lon
  161. REAL :: loninc ! Longitude increment for cylindrical lat/lon
  162. ! also the lon increment for Gaussian grid
  163. REAL :: dlat ! Lat increment for lat/lon grids
  164. REAL :: dlon ! Lon increment for lat/lon grids
  165. REAL :: stdlon ! Longitude parallel to y-axis (-180->180E)
  166. REAL :: truelat1 ! First true latitude (all projections)
  167. REAL :: truelat2 ! Second true lat (LC only)
  168. REAL :: hemi ! 1 for NH, -1 for SH
  169. REAL :: cone ! Cone factor for LC projections
  170. REAL :: polei ! Computed i-location of pole point
  171. REAL :: polej ! Computed j-location of pole point
  172. REAL :: rsw ! Computed radius to SW corner
  173. REAL :: rebydx ! Earth radius divided by dx
  174. REAL :: knowni ! X-location of known lat/lon
  175. REAL :: knownj ! Y-location of known lat/lon
  176. REAL :: re_m ! Radius of spherical earth, meters
  177. REAL :: rho0 ! For Albers equal area
  178. REAL :: nc ! For Albers equal area
  179. REAL :: bigc ! For Albers equal area
  180. LOGICAL :: init ! Flag to indicate if this struct is
  181. ! ready for use
  182. LOGICAL :: wrap ! For Gaussian -- flag to indicate wrapping
  183. ! around globe?
  184. LOGICAL :: comp_ll ! Work in computational lat/lon space for Cassini
  185. REAL, POINTER, DIMENSION(:) :: gauss_lat ! Latitude array for Gaussian grid
  186. END TYPE proj_info
  187. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  188. CONTAINS
  189. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  190. SUBROUTINE map_init(proj)
  191. ! Initializes the map projection structure to missing values
  192. IMPLICIT NONE
  193. TYPE(proj_info), INTENT(INOUT) :: proj
  194. proj%lat1 = -999.9
  195. proj%lon1 = -999.9
  196. proj%lat0 = -999.9
  197. proj%lon0 = -999.9
  198. proj%dx = -999.9
  199. proj%dy = -999.9
  200. proj%latinc = -999.9
  201. proj%loninc = -999.9
  202. proj%stdlon = -999.9
  203. proj%truelat1 = -999.9
  204. proj%truelat2 = -999.9
  205. proj%phi = -999.9
  206. proj%lambda = -999.9
  207. proj%ixdim = -999
  208. proj%jydim = -999
  209. proj%stagger = HH
  210. proj%nlat = 0
  211. proj%nlon = 0
  212. proj%nxmin = 1
  213. proj%nxmax = 43200
  214. proj%hemi = 0.0
  215. proj%cone = -999.9
  216. proj%polei = -999.9
  217. proj%polej = -999.9
  218. proj%rsw = -999.9
  219. proj%knowni = -999.9
  220. proj%knownj = -999.9
  221. proj%re_m = EARTH_RADIUS_M
  222. proj%init = .FALSE.
  223. proj%wrap = .FALSE.
  224. proj%rho0 = 0.
  225. proj%nc = 0.
  226. proj%bigc = 0.
  227. proj%comp_ll = .FALSE.
  228. nullify(proj%gauss_lat)
  229. END SUBROUTINE map_init
  230. SUBROUTINE map_set(proj_code, proj, lat1, lon1, lat0, lon0, knowni, knownj, dx, dy, latinc, &
  231. loninc, stdlon, truelat1, truelat2, nlat, nlon, ixdim, jydim, nxmin, nxmax, &
  232. stagger, phi, lambda, r_earth)
  233. ! Given a partially filled proj_info structure, this routine computes
  234. ! polei, polej, rsw, and cone (if LC projection) to complete the
  235. ! structure. This allows us to eliminate redundant calculations when
  236. ! calling the coordinate conversion routines multiple times for the
  237. ! same map.
  238. ! This will generally be the first routine called when a user wants
  239. ! to be able to use the coordinate conversion routines, and it
  240. ! will call the appropriate subroutines based on the
  241. ! proj%code which indicates which projection type this is.
  242. IMPLICIT NONE
  243. ! Declare arguments
  244. INTEGER, INTENT(IN) :: proj_code
  245. INTEGER, INTENT(IN), OPTIONAL :: nlat
  246. INTEGER, INTENT(IN), OPTIONAL :: nlon
  247. INTEGER, INTENT(IN), OPTIONAL :: ixdim
  248. INTEGER, INTENT(IN), OPTIONAL :: jydim
  249. INTEGER, INTENT(IN), OPTIONAL :: nxmin
  250. INTEGER, INTENT(IN), OPTIONAL :: nxmax
  251. INTEGER, INTENT(IN), OPTIONAL :: stagger
  252. REAL, INTENT(IN), OPTIONAL :: latinc
  253. REAL, INTENT(IN), OPTIONAL :: loninc
  254. REAL, INTENT(IN), OPTIONAL :: lat1
  255. REAL, INTENT(IN), OPTIONAL :: lon1
  256. REAL, INTENT(IN), OPTIONAL :: lat0
  257. REAL, INTENT(IN), OPTIONAL :: lon0
  258. REAL, INTENT(IN), OPTIONAL :: dx
  259. REAL, INTENT(IN), OPTIONAL :: dy
  260. REAL, INTENT(IN), OPTIONAL :: stdlon
  261. REAL, INTENT(IN), OPTIONAL :: truelat1
  262. REAL, INTENT(IN), OPTIONAL :: truelat2
  263. REAL, INTENT(IN), OPTIONAL :: knowni
  264. REAL, INTENT(IN), OPTIONAL :: knownj
  265. REAL, INTENT(IN), OPTIONAL :: phi
  266. REAL, INTENT(IN), OPTIONAL :: lambda
  267. REAL, INTENT(IN), OPTIONAL :: r_earth
  268. TYPE(proj_info), INTENT(OUT) :: proj
  269. INTEGER :: iter
  270. REAL :: dummy_lon1
  271. REAL :: dummy_lon0
  272. REAL :: dummy_stdlon
  273. ! First, verify that mandatory parameters are present for the specified proj_code
  274. IF ( proj_code == PROJ_LC ) THEN
  275. IF ( .NOT.PRESENT(truelat1) .OR. &
  276. .NOT.PRESENT(truelat2) .OR. &
  277. .NOT.PRESENT(lat1) .OR. &
  278. .NOT.PRESENT(lon1) .OR. &
  279. .NOT.PRESENT(knowni) .OR. &
  280. .NOT.PRESENT(knownj) .OR. &
  281. .NOT.PRESENT(stdlon) .OR. &
  282. .NOT.PRESENT(dx) ) THEN
  283. PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
  284. PRINT '(A)', ' truelat1, truelat2, lat1, lon1, knowni, knownj, stdlon, dx'
  285. call mprintf(.true.,ERROR,'MAP_INIT')
  286. END IF
  287. ELSE IF ( proj_code == PROJ_PS ) THEN
  288. IF ( .NOT.PRESENT(truelat1) .OR. &
  289. .NOT.PRESENT(lat1) .OR. &
  290. .NOT.PRESENT(lon1) .OR. &
  291. .NOT.PRESENT(knowni) .OR. &
  292. .NOT.PRESENT(knownj) .OR. &
  293. .NOT.PRESENT(stdlon) .OR. &
  294. .NOT.PRESENT(dx) ) THEN
  295. PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
  296. PRINT '(A)', ' truelat1, lat1, lon1, knonwi, knownj, stdlon, dx'
  297. call mprintf(.true.,ERROR,'MAP_INIT')
  298. END IF
  299. ELSE IF ( proj_code == PROJ_PS_WGS84 ) THEN
  300. IF ( .NOT.PRESENT(truelat1) .OR. &
  301. .NOT.PRESENT(lat1) .OR. &
  302. .NOT.PRESENT(lon1) .OR. &
  303. .NOT.PRESENT(knowni) .OR. &
  304. .NOT.PRESENT(knownj) .OR. &
  305. .NOT.PRESENT(stdlon) .OR. &
  306. .NOT.PRESENT(dx) ) THEN
  307. PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
  308. PRINT '(A)', ' truelat1, lat1, lon1, knonwi, knownj, stdlon, dx'
  309. call mprintf(.true.,ERROR,'MAP_INIT')
  310. END IF
  311. ELSE IF ( proj_code == PROJ_ALBERS_NAD83 ) THEN
  312. IF ( .NOT.PRESENT(truelat1) .OR. &
  313. .NOT.PRESENT(truelat2) .OR. &
  314. .NOT.PRESENT(lat1) .OR. &
  315. .NOT.PRESENT(lon1) .OR. &
  316. .NOT.PRESENT(knowni) .OR. &
  317. .NOT.PRESENT(knownj) .OR. &
  318. .NOT.PRESENT(stdlon) .OR. &
  319. .NOT.PRESENT(dx) ) THEN
  320. PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
  321. PRINT '(A)', ' truelat1, truelat2, lat1, lon1, knonwi, knownj, stdlon, dx'
  322. call mprintf(.true.,ERROR,'MAP_INIT')
  323. END IF
  324. ELSE IF ( proj_code == PROJ_MERC ) THEN
  325. IF ( .NOT.PRESENT(truelat1) .OR. &
  326. .NOT.PRESENT(lat1) .OR. &
  327. .NOT.PRESENT(lon1) .OR. &
  328. .NOT.PRESENT(knowni) .OR. &
  329. .NOT.PRESENT(knownj) .OR. &
  330. .NOT.PRESENT(dx) ) THEN
  331. PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
  332. PRINT '(A)', ' truelat1, lat1, lon1, knowni, knownj, dx'
  333. call mprintf(.true.,ERROR,'MAP_INIT')
  334. END IF
  335. ELSE IF ( proj_code == PROJ_LATLON ) THEN
  336. IF ( .NOT.PRESENT(latinc) .OR. &
  337. .NOT.PRESENT(loninc) .OR. &
  338. .NOT.PRESENT(knowni) .OR. &
  339. .NOT.PRESENT(knownj) .OR. &
  340. .NOT.PRESENT(lat1) .OR. &
  341. .NOT.PRESENT(lon1) ) THEN
  342. PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
  343. PRINT '(A)', ' latinc, loninc, knowni, knownj, lat1, lon1'
  344. call mprintf(.true.,ERROR,'MAP_INIT')
  345. END IF
  346. ELSE IF ( proj_code == PROJ_CYL ) THEN
  347. IF ( .NOT.PRESENT(latinc) .OR. &
  348. .NOT.PRESENT(loninc) .OR. &
  349. .NOT.PRESENT(stdlon) ) THEN
  350. PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
  351. PRINT '(A)', ' latinc, loninc, stdlon'
  352. call mprintf(.true.,ERROR,'MAP_INIT')
  353. END IF
  354. ELSE IF ( proj_code == PROJ_CASSINI ) THEN
  355. IF ( .NOT.PRESENT(latinc) .OR. &
  356. .NOT.PRESENT(loninc) .OR. &
  357. .NOT.PRESENT(lat1) .OR. &
  358. .NOT.PRESENT(lon1) .OR. &
  359. .NOT.PRESENT(lat0) .OR. &
  360. .NOT.PRESENT(lon0) .OR. &
  361. .NOT.PRESENT(knowni) .OR. &
  362. .NOT.PRESENT(knownj) .OR. &
  363. .NOT.PRESENT(stdlon) ) THEN
  364. PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
  365. PRINT '(A)', ' latinc, loninc, lat1, lon1, knowni, knownj, lat0, lon0, stdlon'
  366. call mprintf(.true.,ERROR,'MAP_INIT')
  367. END IF
  368. ELSE IF ( proj_code == PROJ_GAUSS ) THEN
  369. IF ( .NOT.PRESENT(nlat) .OR. &
  370. .NOT.PRESENT(lat1) .OR. &
  371. .NOT.PRESENT(lon1) .OR. &
  372. .NOT.PRESENT(loninc) ) THEN
  373. PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
  374. PRINT '(A)', ' nlat, lat1, lon1, loninc'
  375. call mprintf(.true.,ERROR,'MAP_INIT')
  376. END IF
  377. ELSE IF ( proj_code == PROJ_ROTLL ) THEN
  378. IF ( .NOT.PRESENT(ixdim) .OR. &
  379. .NOT.PRESENT(jydim) .OR. &
  380. .NOT.PRESENT(phi) .OR. &
  381. .NOT.PRESENT(lambda) .OR. &
  382. .NOT.PRESENT(lat1) .OR. &
  383. .NOT.PRESENT(lon1) .OR. &
  384. .NOT.PRESENT(stagger) ) THEN
  385. PRINT '(A,I2)', 'The following are mandatory parameters for projection code : ', proj_code
  386. PRINT '(A)', ' ixdim, jydim, phi, lambda, lat1, lon1, stagger'
  387. call mprintf(.true.,ERROR,'MAP_INIT')
  388. END IF
  389. ELSE
  390. PRINT '(A,I2)', 'Unknown projection code: ', proj_code
  391. call mprintf(.true.,ERROR,'MAP_INIT')
  392. END IF
  393. ! Check for validity of mandatory variables in proj
  394. IF ( PRESENT(lat1) ) THEN
  395. IF ( ABS(lat1) .GT. 90. ) THEN
  396. PRINT '(A)', 'Latitude of origin corner required as follows:'
  397. PRINT '(A)', ' -90N <= lat1 < = 90.N'
  398. call mprintf(.true.,ERROR,'MAP_INIT')
  399. ENDIF
  400. ENDIF
  401. IF ( PRESENT(lon1) ) THEN
  402. dummy_lon1 = lon1
  403. IF ( ABS(dummy_lon1) .GT. 180.) THEN
  404. iter = 0
  405. DO WHILE (ABS(dummy_lon1) > 180. .AND. iter < 10)
  406. IF (dummy_lon1 < -180.) dummy_lon1 = dummy_lon1 + 360.
  407. IF (dummy_lon1 > 180.) dummy_lon1 = dummy_lon1 - 360.
  408. iter = iter + 1
  409. END DO
  410. IF (abs(dummy_lon1) > 180.) THEN
  411. PRINT '(A)', 'Longitude of origin required as follows:'
  412. PRINT '(A)', ' -180E <= lon1 <= 180W'
  413. call mprintf(.true.,ERROR,'MAP_INIT')
  414. ENDIF
  415. ENDIF
  416. ENDIF
  417. IF ( PRESENT(lon0) ) THEN
  418. dummy_lon0 = lon0
  419. IF ( ABS(dummy_lon0) .GT. 180.) THEN
  420. iter = 0
  421. DO WHILE (ABS(dummy_lon0) > 180. .AND. iter < 10)
  422. IF (dummy_lon0 < -180.) dummy_lon0 = dummy_lon0 + 360.
  423. IF (dummy_lon0 > 180.) dummy_lon0 = dummy_lon0 - 360.
  424. iter = iter + 1
  425. END DO
  426. IF (abs(dummy_lon0) > 180.) THEN
  427. PRINT '(A)', 'Longitude of pole required as follows:'
  428. PRINT '(A)', ' -180E <= lon0 <= 180W'
  429. call mprintf(.true.,ERROR,'MAP_INIT')
  430. ENDIF
  431. ENDIF
  432. ENDIF
  433. IF ( PRESENT(dx) ) THEN
  434. IF ((dx .LE. 0.).AND.(proj_code .NE. PROJ_LATLON)) THEN
  435. PRINT '(A)', 'Require grid spacing (dx) in meters be positive!'
  436. call mprintf(.true.,ERROR,'MAP_INIT')
  437. ENDIF
  438. ENDIF
  439. IF ( PRESENT(stdlon) ) THEN
  440. dummy_stdlon = stdlon
  441. IF ((ABS(dummy_stdlon) > 180.).AND.(proj_code /= PROJ_MERC)) THEN
  442. iter = 0
  443. DO WHILE (ABS(dummy_stdlon) > 180. .AND. iter < 10)
  444. IF (dummy_stdlon < -180.) dummy_stdlon = dummy_stdlon + 360.
  445. IF (dummy_stdlon > 180.) dummy_stdlon = dummy_stdlon - 360.
  446. iter = iter + 1
  447. END DO
  448. IF (abs(dummy_stdlon) > 180.) THEN
  449. PRINT '(A)', 'Need orientation longitude (stdlon) as: '
  450. PRINT '(A)', ' -180E <= stdlon <= 180W'
  451. call mprintf(.true.,ERROR,'MAP_INIT')
  452. ENDIF
  453. ENDIF
  454. ENDIF
  455. IF ( PRESENT(truelat1) ) THEN
  456. IF (ABS(truelat1).GT.90.) THEN
  457. PRINT '(A)', 'Set true latitude 1 for all projections!'
  458. call mprintf(.true.,ERROR,'MAP_INIT')
  459. ENDIF
  460. ENDIF
  461. CALL map_init(proj)
  462. proj%code = proj_code
  463. IF ( PRESENT(lat1) ) proj%lat1 = lat1
  464. IF ( PRESENT(lon1) ) proj%lon1 = dummy_lon1
  465. IF ( PRESENT(lat0) ) proj%lat0 = lat0
  466. IF ( PRESENT(lon0) ) proj%lon0 = dummy_lon0
  467. IF ( PRESENT(latinc) ) proj%latinc = latinc
  468. IF ( PRESENT(loninc) ) proj%loninc = loninc
  469. IF ( PRESENT(knowni) ) proj%knowni = knowni
  470. IF ( PRESENT(knownj) ) proj%knownj = knownj
  471. IF ( PRESENT(nxmin) ) proj%nxmin = nxmin
  472. IF ( PRESENT(nxmax) ) proj%nxmax = nxmax
  473. IF ( PRESENT(dx) ) proj%dx = dx
  474. IF ( PRESENT(dy) ) THEN
  475. proj%dy = dy
  476. ELSE IF ( PRESENT(dx) ) THEN
  477. proj%dy = dx
  478. END IF
  479. IF ( PRESENT(stdlon) ) proj%stdlon = dummy_stdlon
  480. IF ( PRESENT(truelat1) ) proj%truelat1 = truelat1
  481. IF ( PRESENT(truelat2) ) proj%truelat2 = truelat2
  482. IF ( PRESENT(nlat) ) proj%nlat = nlat
  483. IF ( PRESENT(nlon) ) proj%nlon = nlon
  484. IF ( PRESENT(ixdim) ) proj%ixdim = ixdim
  485. IF ( PRESENT(jydim) ) proj%jydim = jydim
  486. IF ( PRESENT(stagger) ) proj%stagger = stagger
  487. IF ( PRESENT(phi) ) proj%phi = phi
  488. IF ( PRESENT(lambda) ) proj%lambda = lambda
  489. IF ( PRESENT(r_earth) ) proj%re_m = r_earth
  490. IF ( PRESENT(dx) ) THEN
  491. IF ( (proj_code == PROJ_LC) .OR. (proj_code == PROJ_PS) .OR. &
  492. (proj_code == PROJ_PS_WGS84) .OR. (proj_code == PROJ_ALBERS_NAD83) .OR. &
  493. (proj_code == PROJ_MERC) ) THEN
  494. proj%dx = dx
  495. IF (truelat1 .LT. 0.) THEN
  496. proj%hemi = -1.0
  497. ELSE
  498. proj%hemi = 1.0
  499. ENDIF
  500. proj%rebydx = proj%re_m / dx
  501. ENDIF
  502. ENDIF
  503. pick_proj: SELECT CASE(proj%code)
  504. CASE(PROJ_PS)
  505. CALL set_ps(proj)
  506. CASE(PROJ_PS_WGS84)
  507. CALL set_ps_wgs84(proj)
  508. CASE(PROJ_ALBERS_NAD83)
  509. CALL set_albers_nad83(proj)
  510. CASE(PROJ_LC)
  511. IF (ABS(proj%truelat2) .GT. 90.) THEN
  512. proj%truelat2=proj%truelat1
  513. ENDIF
  514. CALL set_lc(proj)
  515. CASE (PROJ_MERC)
  516. CALL set_merc(proj)
  517. CASE (PROJ_LATLON)
  518. CASE (PROJ_GAUSS)
  519. CALL set_gauss(proj)
  520. CASE (PROJ_CYL)
  521. CALL set_cyl(proj)
  522. CASE (PROJ_CASSINI)
  523. CALL set_cassini(proj)
  524. CASE (PROJ_ROTLL)
  525. END SELECT pick_proj
  526. proj%init = .TRUE.
  527. RETURN
  528. END SUBROUTINE map_set
  529. SUBROUTINE latlon_to_ij(proj, lat, lon, i, j)
  530. ! Converts input lat/lon values to the cartesian (i,j) value
  531. ! for the given projection.
  532. IMPLICIT NONE
  533. TYPE(proj_info), INTENT(IN) :: proj
  534. REAL, INTENT(IN) :: lat
  535. REAL, INTENT(IN) :: lon
  536. REAL, INTENT(OUT) :: i
  537. REAL, INTENT(OUT) :: j
  538. IF (.NOT.proj%init) THEN
  539. PRINT '(A)', 'You have not called map_set for this projection!'
  540. call mprintf(.true.,ERROR,'LATLON_TO_IJ')
  541. ENDIF
  542. SELECT CASE(proj%code)
  543. CASE(PROJ_LATLON)
  544. CALL llij_latlon(lat,lon,proj,i,j)
  545. CASE(PROJ_MERC)
  546. CALL llij_merc(lat,lon,proj,i,j)
  547. CASE(PROJ_PS)
  548. CALL llij_ps(lat,lon,proj,i,j)
  549. CASE(PROJ_PS_WGS84)
  550. CALL llij_ps_wgs84(lat,lon,proj,i,j)
  551. CASE(PROJ_ALBERS_NAD83)
  552. CALL llij_albers_nad83(lat,lon,proj,i,j)
  553. CASE(PROJ_LC)
  554. CALL llij_lc(lat,lon,proj,i,j)
  555. CASE(PROJ_GAUSS)
  556. CALL llij_gauss(lat,lon,proj,i,j)
  557. CASE(PROJ_CYL)
  558. CALL llij_cyl(lat,lon,proj,i,j)
  559. CASE(PROJ_CASSINI)
  560. CALL llij_cassini(lat,lon,proj,i,j)
  561. CASE(PROJ_ROTLL)
  562. CALL llij_rotlatlon(lat,lon,proj,i,j)
  563. CASE DEFAULT
  564. PRINT '(A,I2)', 'Unrecognized map projection code: ', proj%code
  565. call mprintf(.true.,ERROR,'LATLON_TO_IJ')
  566. END SELECT
  567. RETURN
  568. END SUBROUTINE latlon_to_ij
  569. SUBROUTINE ij_to_latlon(proj, i, j, lat, lon)
  570. ! Computes geographical latitude and longitude for a given (i,j) point
  571. ! in a grid with a projection of proj
  572. IMPLICIT NONE
  573. TYPE(proj_info),INTENT(IN) :: proj
  574. REAL, INTENT(IN) :: i
  575. REAL, INTENT(IN) :: j
  576. REAL, INTENT(OUT) :: lat
  577. REAL, INTENT(OUT) :: lon
  578. IF (.NOT.proj%init) THEN
  579. PRINT '(A)', 'You have not called map_set for this projection!'
  580. call mprintf(.true.,ERROR,'IJ_TO_LATLON')
  581. ENDIF
  582. SELECT CASE (proj%code)
  583. CASE (PROJ_LATLON)
  584. CALL ijll_latlon(i, j, proj, lat, lon)
  585. CASE (PROJ_MERC)
  586. CALL ijll_merc(i, j, proj, lat, lon)
  587. CASE (PROJ_PS)
  588. CALL ijll_ps(i, j, proj, lat, lon)
  589. CASE (PROJ_PS_WGS84)
  590. CALL ijll_ps_wgs84(i, j, proj, lat, lon)
  591. CASE (PROJ_ALBERS_NAD83)
  592. CALL ijll_albers_nad83(i, j, proj, lat, lon)
  593. CASE (PROJ_LC)
  594. CALL ijll_lc(i, j, proj, lat, lon)
  595. CASE (PROJ_CYL)
  596. CALL ijll_cyl(i, j, proj, lat, lon)
  597. CASE (PROJ_CASSINI)
  598. CALL ijll_cassini(i, j, proj, lat, lon)
  599. CASE (PROJ_ROTLL)
  600. CALL ijll_rotlatlon(i, j, proj, lat, lon)
  601. CASE DEFAULT
  602. PRINT '(A,I2)', 'Unrecognized map projection code: ', proj%code
  603. call mprintf(.true.,ERROR,'IJ_TO_LATLON')
  604. END SELECT
  605. RETURN
  606. END SUBROUTINE ij_to_latlon
  607. SUBROUTINE set_ps(proj)
  608. ! Initializes a polar-stereographic map projection from the partially
  609. ! filled proj structure. This routine computes the radius to the
  610. ! southwest corner and computes the i/j location of the pole for use
  611. ! in llij_ps and ijll_ps.
  612. IMPLICIT NONE
  613. ! Declare args
  614. TYPE(proj_info), INTENT(INOUT) :: proj
  615. ! Local vars
  616. REAL :: ala1
  617. REAL :: alo1
  618. REAL :: reflon
  619. REAL :: scale_top
  620. ! Executable code
  621. reflon = proj%stdlon + 90.
  622. ! Compute numerator term of map scale factor
  623. scale_top = 1. + proj%hemi * SIN(proj%truelat1 * rad_per_deg)
  624. ! Compute radius to lower-left (SW) corner
  625. ala1 = proj%lat1 * rad_per_deg
  626. proj%rsw = proj%rebydx*COS(ala1)*scale_top/(1.+proj%hemi*SIN(ala1))
  627. ! Find the pole point
  628. alo1 = (proj%lon1 - reflon) * rad_per_deg
  629. proj%polei = proj%knowni - proj%rsw * COS(alo1)
  630. proj%polej = proj%knownj - proj%hemi * proj%rsw * SIN(alo1)
  631. RETURN
  632. END SUBROUTINE set_ps
  633. SUBROUTINE llij_ps(lat,lon,proj,i,j)
  634. ! Given latitude (-90 to 90), longitude (-180 to 180), and the
  635. ! standard polar-stereographic projection information via the
  636. ! public proj structure, this routine returns the i/j indices which
  637. ! if within the domain range from 1->nx and 1->ny, respectively.
  638. IMPLICIT NONE
  639. ! Delcare input arguments
  640. REAL, INTENT(IN) :: lat
  641. REAL, INTENT(IN) :: lon
  642. TYPE(proj_info),INTENT(IN) :: proj
  643. ! Declare output arguments
  644. REAL, INTENT(OUT) :: i !(x-index)
  645. REAL, INTENT(OUT) :: j !(y-index)
  646. ! Declare local variables
  647. REAL :: reflon
  648. REAL :: scale_top
  649. REAL :: ala
  650. REAL :: alo
  651. REAL :: rm
  652. ! BEGIN CODE
  653. reflon = proj%stdlon + 90.
  654. ! Compute numerator term of map scale factor
  655. scale_top = 1. + proj%hemi * SIN(proj%truelat1 * rad_per_deg)
  656. ! Find radius to desired point
  657. ala = lat * rad_per_deg
  658. rm = proj%rebydx * COS(ala) * scale_top/(1. + proj%hemi *SIN(ala))
  659. alo = (lon - reflon) * rad_per_deg
  660. i = proj%polei + rm * COS(alo)
  661. j = proj%polej + proj%hemi * rm * SIN(alo)
  662. RETURN
  663. END SUBROUTINE llij_ps
  664. SUBROUTINE ijll_ps(i, j, proj, lat, lon)
  665. ! This is the inverse subroutine of llij_ps. It returns the
  666. ! latitude and longitude of an i/j point given the projection info
  667. ! structure.
  668. IMPLICIT NONE
  669. ! Declare input arguments
  670. REAL, INTENT(IN) :: i ! Column
  671. REAL, INTENT(IN) :: j ! Row
  672. TYPE (proj_info), INTENT(IN) :: proj
  673. ! Declare output arguments
  674. REAL, INTENT(OUT) :: lat ! -90 -> 90 north
  675. REAL, INTENT(OUT) :: lon ! -180 -> 180 East
  676. ! Local variables
  677. REAL :: reflon
  678. REAL :: scale_top
  679. REAL :: xx,yy
  680. REAL :: gi2, r2
  681. REAL :: arccos
  682. ! Begin Code
  683. ! Compute the reference longitude by rotating 90 degrees to the east
  684. ! to find the longitude line parallel to the positive x-axis.
  685. reflon = proj%stdlon + 90.
  686. ! Compute numerator term of map scale factor
  687. scale_top = 1. + proj%hemi * SIN(proj%truelat1 * rad_per_deg)
  688. ! Compute radius to point of interest
  689. xx = i - proj%polei
  690. yy = (j - proj%polej) * proj%hemi
  691. r2 = xx**2 + yy**2
  692. ! Now the magic code
  693. IF (r2 .EQ. 0.) THEN
  694. lat = proj%hemi * 90.
  695. lon = reflon
  696. ELSE
  697. gi2 = (proj%rebydx * scale_top)**2.
  698. lat = deg_per_rad * proj%hemi * ASIN((gi2-r2)/(gi2+r2))
  699. arccos = ACOS(xx/SQRT(r2))
  700. IF (yy .GT. 0) THEN
  701. lon = reflon + deg_per_rad * arccos
  702. ELSE
  703. lon = reflon - deg_per_rad * arccos
  704. ENDIF
  705. ENDIF
  706. ! Convert to a -180 -> 180 East convention
  707. IF (lon .GT. 180.) lon = lon - 360.
  708. IF (lon .LT. -180.) lon = lon + 360.
  709. RETURN
  710. END SUBROUTINE ijll_ps
  711. SUBROUTINE set_ps_wgs84(proj)
  712. ! Initializes a polar-stereographic map projection (WGS84 ellipsoid)
  713. ! from the partially filled proj structure. This routine computes the
  714. ! radius to the southwest corner and computes the i/j location of the
  715. ! pole for use in llij_ps and ijll_ps.
  716. IMPLICIT NONE
  717. ! Arguments
  718. TYPE(proj_info), INTENT(INOUT) :: proj
  719. ! Local variables
  720. real :: h, mc, tc, t, rho
  721. h = proj%hemi
  722. mc = cos(h*proj%truelat1*rad_per_deg)/sqrt(1.0-(E_WGS84*sin(h*proj%truelat1*rad_per_deg))**2.0)
  723. tc = sqrt(((1.0-sin(h*proj%truelat1*rad_per_deg))/(1.0+sin(h*proj%truelat1*rad_per_deg)))* &
  724. (((1.0+E_WGS84*sin(h*proj%truelat1*rad_per_deg))/(1.0-E_WGS84*sin(h*proj%truelat1*rad_per_deg)))**E_WGS84 ))
  725. ! Find the i/j location of reference lat/lon with respect to the pole of the projection
  726. t = sqrt(((1.0-sin(h*proj%lat1*rad_per_deg))/(1.0+sin(h*proj%lat1*rad_per_deg)))* &
  727. (((1.0+E_WGS84*sin(h*proj%lat1*rad_per_deg))/(1.0-E_WGS84*sin(h*proj%lat1*rad_per_deg)) )**E_WGS84 ) )
  728. rho = h * (A_WGS84 / proj%dx) * mc * t / tc
  729. proj%polei = rho * sin((h*proj%lon1 - h*proj%stdlon)*rad_per_deg)
  730. proj%polej = -rho * cos((h*proj%lon1 - h*proj%stdlon)*rad_per_deg)
  731. RETURN
  732. END SUBROUTINE set_ps_wgs84
  733. SUBROUTINE llij_ps_wgs84(lat,lon,proj,i,j)
  734. ! Given latitude (-90 to 90), longitude (-180 to 180), and the
  735. ! standard polar-stereographic projection information via the
  736. ! public proj structure, this routine returns the i/j indices which
  737. ! if within the domain range from 1->nx and 1->ny, respectively.
  738. IMPLICIT NONE
  739. ! Arguments
  740. REAL, INTENT(IN) :: lat
  741. REAL, INTENT(IN) :: lon
  742. REAL, INTENT(OUT) :: i !(x-index)
  743. REAL, INTENT(OUT) :: j !(y-index)
  744. TYPE(proj_info),INTENT(IN) :: proj
  745. ! Local variables
  746. real :: h, mc, tc, t, rho
  747. h = proj%hemi
  748. mc = cos(h*proj%truelat1*rad_per_deg)/sqrt(1.0-(E_WGS84*sin(h*proj%truelat1*rad_per_deg))**2.0)
  749. tc = sqrt(((1.0-sin(h*proj%truelat1*rad_per_deg))/(1.0+sin(h*proj%truelat1*rad_per_deg)))* &
  750. (((1.0+E_WGS84*sin(h*proj%truelat1*rad_per_deg))/(1.0-E_WGS84*sin(h*proj%truelat1*rad_per_deg)))**E_WGS84 ))
  751. t = sqrt(((1.0-sin(h*lat*rad_per_deg))/(1.0+sin(h*lat*rad_per_deg))) * &
  752. (((1.0+E_WGS84*sin(h*lat*rad_per_deg))/(1.0-E_WGS84*sin(h*lat*rad_per_deg)))**E_WGS84))
  753. ! Find the x/y location of the requested lat/lon with respect to the pole of the projection
  754. rho = (A_WGS84 / proj%dx) * mc * t / tc
  755. i = h * rho * sin((h*lon - h*proj%stdlon)*rad_per_deg)
  756. j = h *(-rho)* cos((h*lon - h*proj%stdlon)*rad_per_deg)
  757. ! Get i/j relative to reference i/j
  758. i = proj%knowni + (i - proj%polei)
  759. j = proj%knownj + (j - proj%polej)
  760. RETURN
  761. END SUBROUTINE llij_ps_wgs84
  762. SUBROUTINE ijll_ps_wgs84(i, j, proj, lat, lon)
  763. ! This is the inverse subroutine of llij_ps. It returns the
  764. ! latitude and longitude of an i/j point given the projection info
  765. ! structure.
  766. implicit none
  767. ! Arguments
  768. REAL, INTENT(IN) :: i ! Column
  769. REAL, INTENT(IN) :: j ! Row
  770. REAL, INTENT(OUT) :: lat ! -90 -> 90 north
  771. REAL, INTENT(OUT) :: lon ! -180 -> 180 East
  772. TYPE (proj_info), INTENT(IN) :: proj
  773. ! Local variables
  774. real :: h, mc, tc, t, rho, x, y
  775. real :: chi, a, b, c, d
  776. h = proj%hemi
  777. x = (i - proj%knowni + proj%polei)
  778. y = (j - proj%knownj + proj%polej)
  779. mc = cos(h*proj%truelat1*rad_per_deg)/sqrt(1.0-(E_WGS84*sin(h*proj%truelat1*rad_per_deg))**2.0)
  780. tc = sqrt(((1.0-sin(h*proj%truelat1*rad_per_deg))/(1.0+sin(h*proj%truelat1*rad_per_deg))) * &
  781. (((1.0+E_WGS84*sin(h*proj%truelat1*rad_per_deg))/(1.0-E_WGS84*sin(h*proj%truelat1*rad_per_deg)))**E_WGS84 ))
  782. rho = sqrt((x*proj%dx)**2.0 + (y*proj%dx)**2.0)
  783. t = rho * tc / (A_WGS84 * mc)
  784. lon = h*proj%stdlon + h*atan2(h*x,h*(-y))
  785. chi = PI/2.0-2.0*atan(t)
  786. a = 1./2.*E_WGS84**2. + 5./24.*E_WGS84**4. + 1./40.*E_WGS84**6. + 73./2016.*E_WGS84**8.
  787. b = 7./24.*E_WGS84**4. + 29./120.*E_WGS84**6. + 54113./40320.*E_WGS84**8.
  788. c = 7./30.*E_WGS84**6. + 81./280.*E_WGS84**8.
  789. d = 4279./20160.*E_WGS84**8.
  790. lat = chi + sin(2.*chi)*(a + cos(2.*chi)*(b + cos(2.*chi)*(c + d*cos(2.*chi))))
  791. lat = h * lat
  792. lat = lat*deg_per_rad
  793. lon = lon*deg_per_rad
  794. RETURN
  795. END SUBROUTINE ijll_ps_wgs84
  796. SUBROUTINE set_albers_nad83(proj)
  797. ! Initializes an Albers equal area map projection (NAD83 ellipsoid)
  798. ! from the partially filled proj structure. This routine computes the
  799. ! radius to the southwest corner and computes the i/j location of the
  800. ! pole for use in llij_albers_nad83 and ijll_albers_nad83.
  801. IMPLICIT NONE
  802. ! Arguments
  803. TYPE(proj_info), INTENT(INOUT) :: proj
  804. ! Local variables
  805. real :: h, m1, m2, q1, q2, theta, q, sinphi
  806. h = proj%hemi
  807. m1 = cos(h*proj%truelat1*rad_per_deg)/sqrt(1.0-(E_NAD83*sin(h*proj%truelat1*rad_per_deg))**2.0)
  808. m2 = cos(h*proj%truelat2*rad_per_deg)/sqrt(1.0-(E_NAD83*sin(h*proj%truelat2*rad_per_deg))**2.0)
  809. sinphi = sin(proj%truelat1*rad_per_deg)
  810. q1 = (1.0-E_NAD83**2.0) * &
  811. ((sinphi/(1.0-(E_NAD83*sinphi)**2.0)) - 1.0/(2.0*E_NAD83) * log((1.0-E_NAD83*sinphi)/(1.0+E_NAD83*sinphi)))
  812. sinphi = sin(proj%truelat2*rad_per_deg)
  813. q2 = (1.0-E_NAD83**2.0) * &
  814. ((sinphi/(1.0-(E_NAD83*sinphi)**2.0)) - 1.0/(2.0*E_NAD83) * log((1.0-E_NAD83*sinphi)/(1.0+E_NAD83*sinphi)))
  815. if (proj%truelat1 == proj%truelat2) then
  816. proj%nc = sin(proj%truelat1*rad_per_deg)
  817. else
  818. proj%nc = (m1**2.0 - m2**2.0) / (q2 - q1)
  819. end if
  820. proj%bigc = m1**2.0 + proj%nc*q1
  821. ! Find the i/j location of reference lat/lon with respect to the pole of the projection
  822. sinphi = sin(proj%lat1*rad_per_deg)
  823. q = (1.0-E_NAD83**2.0) * &
  824. ((sinphi/(1.0-(E_NAD83*sinphi)**2.0)) - 1.0/(2.0*E_NAD83) * log((1.0-E_NAD83*sinphi)/(1.0+E_NAD83*sinphi)))
  825. proj%rho0 = h * (A_NAD83 / proj%dx) * sqrt(proj%bigc - proj%nc * q) / proj%nc
  826. theta = proj%nc*(proj%lon1 - proj%stdlon)*rad_per_deg
  827. proj%polei = proj%rho0 * sin(h*theta)
  828. proj%polej = proj%rho0 - proj%rho0 * cos(h*theta)
  829. RETURN
  830. END SUBROUTINE set_albers_nad83
  831. SUBROUTINE llij_albers_nad83(lat,lon,proj,i,j)
  832. ! Given latitude (-90 to 90), longitude (-180 to 180), and the
  833. ! standard projection information via the
  834. ! public proj structure, this routine returns the i/j indices which
  835. ! if within the domain range from 1->nx and 1->ny, respectively.
  836. IMPLICIT NONE
  837. ! Arguments
  838. REAL, INTENT(IN) :: lat
  839. REAL, INTENT(IN) :: lon
  840. REAL, INTENT(OUT) :: i !(x-index)
  841. REAL, INTENT(OUT) :: j !(y-index)
  842. TYPE(proj_info),INTENT(IN) :: proj
  843. ! Local variables
  844. real :: h, q, rho, theta, sinphi
  845. h = proj%hemi
  846. sinphi = sin(h*lat*rad_per_deg)
  847. ! Find the x/y location of the requested lat/lon with respect to the pole of the projection
  848. q = (1.0-E_NAD83**2.0) * &
  849. ((sinphi/(1.0-(E_NAD83*sinphi)**2.0)) - 1.0/(2.0*E_NAD83) * log((1.0-E_NAD83*sinphi)/(1.0+E_NAD83*sinphi)))
  850. rho = h * (A_NAD83 / proj%dx) * sqrt(proj%bigc - proj%nc * q) / proj%nc
  851. theta = proj%nc * (h*lon - h*proj%stdlon)*rad_per_deg
  852. i = h*rho*sin(theta)
  853. j = h*proj%rho0 - h*rho*cos(theta)
  854. ! Get i/j relative to reference i/j
  855. i = proj%knowni + (i - proj%polei)
  856. j = proj%knownj + (j - proj%polej)
  857. RETURN
  858. END SUBROUTINE llij_albers_nad83
  859. SUBROUTINE ijll_albers_nad83(i, j, proj, lat, lon)
  860. ! This is the inverse subroutine of llij_albers_nad83. It returns the
  861. ! latitude and longitude of an i/j point given the projection info
  862. ! structure.
  863. implicit none
  864. ! Arguments
  865. REAL, INTENT(IN) :: i ! Column
  866. REAL, INTENT(IN) :: j ! Row
  867. REAL, INTENT(OUT) :: lat ! -90 -> 90 north
  868. REAL, INTENT(OUT) :: lon ! -180 -> 180 East
  869. TYPE (proj_info), INTENT(IN) :: proj
  870. ! Local variables
  871. real :: h, q, rho, theta, beta, x, y
  872. real :: a, b, c
  873. h = proj%hemi
  874. x = (i - proj%knowni + proj%polei)
  875. y = (j - proj%knownj + proj%polej)
  876. rho = sqrt(x**2.0 + (proj%rho0 - y)**2.0)
  877. theta = atan2(x, proj%rho0-y)
  878. q = (proj%bigc - (rho*proj%nc*proj%dx/A_NAD83)**2.0) / proj%nc
  879. beta = asin(q/(1.0 - log((1.0-E_NAD83)/(1.0+E_NAD83))*(1.0-E_NAD83**2.0)/(2.0*E_NAD83)))
  880. a = 1./3.*E_NAD83**2. + 31./180.*E_NAD83**4. + 517./5040.*E_NAD83**6.
  881. b = 23./360.*E_NAD83**4. + 251./3780.*E_NAD83**6.
  882. c = 761./45360.*E_NAD83**6.
  883. lat = beta + a*sin(2.*beta) + b*sin(4.*beta) + c*sin(6.*beta)
  884. lat = h*lat*deg_per_rad
  885. lon = proj%stdlon + theta*deg_per_rad/proj%nc
  886. RETURN
  887. END SUBROUTINE ijll_albers_nad83
  888. SUBROUTINE set_lc(proj)
  889. ! Initialize the remaining items in the proj structure for a
  890. ! lambert conformal grid.
  891. IMPLICIT NONE
  892. TYPE(proj_info), INTENT(INOUT) :: proj
  893. REAL :: arg
  894. REAL :: deltalon1
  895. REAL :: tl1r
  896. REAL :: ctl1r
  897. ! Compute cone factor
  898. CALL lc_cone(proj%truelat1, proj%truelat2, proj%cone)
  899. ! Compute longitude differences and ensure we stay out of the
  900. ! forbidden "cut zone"
  901. deltalon1 = proj%lon1 - proj%stdlon
  902. IF (deltalon1 .GT. +180.) deltalon1 = deltalon1 - 360.
  903. IF (deltalon1 .LT. -180.) deltalon1 = deltalon1 + 360.
  904. ! Convert truelat1 to radian and compute COS for later use
  905. tl1r = proj%truelat1 * rad_per_deg
  906. ctl1r = COS(tl1r)
  907. ! Compute the radius to our known lower-left (SW) corner
  908. proj%rsw = proj%rebydx * ctl1r/proj%cone * &
  909. (TAN((90.*proj%hemi-proj%lat1)*rad_per_deg/2.) / &
  910. TAN((90.*proj%hemi-proj%truelat1)*rad_per_deg/2.))**proj%cone
  911. ! Find pole point
  912. arg = proj%cone*(deltalon1*rad_per_deg)
  913. proj%polei = proj%hemi*proj%knowni - proj%hemi * proj%rsw * SIN(arg)
  914. proj%polej = proj%hemi*proj%knownj + proj%rsw * COS(arg)
  915. RETURN
  916. END SUBROUTINE set_lc
  917. SUBROUTINE lc_cone(truelat1, truelat2, cone)
  918. ! Subroutine to compute the cone factor of a Lambert Conformal projection
  919. IMPLICIT NONE
  920. ! Input Args
  921. REAL, INTENT(IN) :: truelat1 ! (-90 -> 90 degrees N)
  922. REAL, INTENT(IN) :: truelat2 ! " " " " "
  923. ! Output Args
  924. REAL, INTENT(OUT) :: cone
  925. ! Locals
  926. ! BEGIN CODE
  927. ! First, see if this is a secant or tangent projection. For tangent
  928. ! projections, truelat1 = truelat2 and the cone is tangent to the
  929. ! Earth's surface at this latitude. For secant projections, the cone
  930. ! intersects the Earth's surface at each of the distinctly different
  931. ! latitudes
  932. IF (ABS(truelat1-truelat2) .GT. 0.1) THEN
  933. cone = ALOG10(COS(truelat1*rad_per_deg)) - &
  934. ALOG10(COS(truelat2*rad_per_deg))
  935. cone = cone /(ALOG10(TAN((45.0 - ABS(truelat1)/2.0) * rad_per_deg)) - &
  936. ALOG10(TAN((45.0 - ABS(truelat2)/2.0) * rad_per_deg)))
  937. ELSE
  938. cone = SIN(ABS(truelat1)*rad_per_deg )
  939. ENDIF
  940. RETURN
  941. END SUBROUTINE lc_cone
  942. SUBROUTINE ijll_lc( i, j, proj, lat, lon)
  943. ! Subroutine to convert from the (i,j) cartesian coordinate to the
  944. ! geographical latitude and longitude for a Lambert Conformal projection.
  945. ! History:
  946. ! 25 Jul 01: Corrected by B. Shaw, NOAA/FSL
  947. !
  948. IMPLICIT NONE
  949. ! Input Args
  950. REAL, INTENT(IN) :: i ! Cartesian X coordinate
  951. REAL, INTENT(IN) :: j ! Cartesian Y coordinate
  952. TYPE(proj_info),INTENT(IN) :: proj ! Projection info structure
  953. ! Output Args
  954. REAL, INTENT(OUT) :: lat ! Latitude (-90->90 deg N)
  955. REAL, INTENT(OUT) :: lon ! Longitude (-180->180 E)
  956. ! Locals
  957. REAL :: inew
  958. REAL :: jnew
  959. REAL :: r
  960. REAL :: chi,chi1,chi2
  961. REAL :: r2
  962. REAL :: xx
  963. REAL :: yy
  964. ! BEGIN CODE
  965. chi1 = (90. - proj%hemi*proj%truelat1)*rad_per_deg
  966. chi2 = (90. - proj%hemi*proj%truelat2)*rad_per_deg
  967. ! See if we are in the southern hemispere and flip the indices
  968. ! if we are.
  969. inew = proj%hemi * i
  970. jnew = proj%hemi * j
  971. ! Compute radius**2 to i/j location
  972. xx = inew - proj%polei
  973. yy = proj%polej - jnew
  974. r2 = (xx*xx + yy*yy)
  975. r = SQRT(r2)/proj%rebydx
  976. ! Convert to lat/lon
  977. IF (r2 .EQ. 0.) THEN
  978. lat = proj%hemi * 90.
  979. lon = proj%stdlon
  980. ELSE
  981. ! Longitude
  982. lon = proj%stdlon + deg_per_rad * ATAN2(proj%hemi*xx,yy)/proj%cone
  983. lon = AMOD(lon+360., 360.)
  984. ! Latitude. Latitude determined by solving an equation adapted
  985. ! from:
  986. ! Maling, D.H., 1973: Coordinate Systems and Map Projections
  987. ! Equations #20 in Appendix I.
  988. IF (chi1 .EQ. chi2) THEN
  989. chi = 2.0*ATAN( ( r/TAN(chi1) )**(1./proj%cone) * TAN(chi1*0.5) )
  990. ELSE
  991. chi = 2.0*ATAN( (r*proj%cone/SIN(chi1))**(1./proj%cone) * TAN(chi1*0.5))
  992. ENDIF
  993. lat = (90.0-chi*deg_per_rad)*proj%hemi
  994. ENDIF
  995. IF (lon .GT. +180.) lon = lon - 360.
  996. IF (lon .LT. -180.) lon = lon + 360.
  997. RETURN
  998. END SUBROUTINE ijll_lc
  999. SUBROUTINE llij_lc( lat, lon, proj, i, j)
  1000. ! Subroutine to compute the geographical latitude and longitude values
  1001. ! to the cartesian x/y on a Lambert Conformal projection.
  1002. IMPLICIT NONE
  1003. ! Input Args
  1004. REAL, INTENT(IN) :: lat ! Latitude (-90->90 deg N)
  1005. REAL, INTENT(IN) :: lon ! Longitude (-180->180 E)
  1006. TYPE(proj_info),INTENT(IN) :: proj ! Projection info structure
  1007. ! Output Args
  1008. REAL, INTENT(OUT) :: i ! Cartesian X coordinate
  1009. REAL, INTENT(OUT) :: j ! Cartesian Y coordinate
  1010. ! Locals
  1011. REAL :: arg
  1012. REAL :: deltalon
  1013. REAL :: tl1r
  1014. REAL :: rm
  1015. REAL :: ctl1r
  1016. ! BEGIN CODE
  1017. ! Compute deltalon between known longitude and standard lon and ensure
  1018. ! it is not in the cut zone
  1019. deltalon = lon - proj%stdlon
  1020. IF (deltalon .GT. +180.) deltalon = deltalon - 360.
  1021. IF (deltalon .LT. -180.) deltalon = deltalon + 360.
  1022. ! Convert truelat1 to radian and compute COS for later use
  1023. tl1r = proj%truelat1 * rad_per_deg
  1024. ctl1r = COS(tl1r)
  1025. ! Radius to desired point
  1026. rm = proj%rebydx * ctl1r/proj%cone * &
  1027. (TAN((90.*proj%hemi-lat)*rad_per_deg/2.) / &
  1028. TAN((90.*proj%hemi-proj%truelat1)*rad_per_deg/2.))**proj%cone
  1029. arg = proj%cone*(deltalon*rad_per_deg)
  1030. i = proj%polei + proj%hemi * rm * SIN(arg)
  1031. j = proj%polej - rm * COS(arg)
  1032. ! Finally, if we are in the southern hemisphere, flip the i/j
  1033. ! values to a coordinate system where (1,1) is the SW corner
  1034. ! (what we assume) which is different than the original NCEP
  1035. ! algorithms which used the NE corner as the origin in the
  1036. ! southern hemisphere (left-hand vs. right-hand coordinate?)
  1037. i = proj%hemi * i
  1038. j = proj%hemi * j
  1039. RETURN
  1040. END SUBROUTINE llij_lc
  1041. SUBROUTINE set_merc(proj)
  1042. ! Sets up the remaining basic elements for the mercator projection
  1043. IMPLICIT NONE
  1044. TYPE(proj_info), INTENT(INOUT) :: proj
  1045. REAL :: clain
  1046. ! Preliminary variables
  1047. clain = COS(rad_per_deg*proj%truelat1)
  1048. proj%dlon = proj%dx / (proj%re_m * clain)
  1049. ! Compute distance from equator to origin, and store in the
  1050. ! proj%rsw tag.
  1051. proj%rsw = 0.
  1052. IF (proj%lat1 .NE. 0.) THEN
  1053. proj%rsw = (ALOG(TAN(0.5*((proj%lat1+90.)*rad_per_deg))))/proj%dlon
  1054. ENDIF
  1055. RETURN
  1056. END SUBROUTINE set_merc
  1057. SUBROUTINE llij_merc(lat, lon, proj, i, j)
  1058. ! Compute i/j coordinate from lat lon for mercator projection
  1059. IMPLICIT NONE
  1060. REAL, INTENT(IN) :: lat
  1061. REAL, INTENT(IN) :: lon
  1062. TYPE(proj_info),INTENT(IN) :: proj
  1063. REAL,INTENT(OUT) :: i
  1064. REAL,INTENT(OUT) :: j
  1065. REAL :: deltalon
  1066. deltalon = lon - proj%lon1
  1067. IF (deltalon .LT. -180.) deltalon = deltalon + 360.
  1068. IF (deltalon .GT. 180.) deltalon = deltalon - 360.
  1069. i = proj%knowni + (deltalon/(proj%dlon*deg_per_rad))
  1070. j = proj%knownj + (ALOG(TAN(0.5*((lat + 90.) * rad_per_deg)))) / &
  1071. proj%dlon - proj%rsw
  1072. RETURN
  1073. END SUBROUTINE llij_merc
  1074. SUBROUTINE ijll_merc(i, j, proj, lat, lon)
  1075. ! Compute the lat/lon from i/j for mercator projection
  1076. IMPLICIT NONE
  1077. REAL,INTENT(IN) :: i
  1078. REAL,INTENT(IN) :: j
  1079. TYPE(proj_info),INTENT(IN) :: proj
  1080. REAL, INTENT(OUT) :: lat
  1081. REAL, INTENT(OUT) :: lon
  1082. lat = 2.0*ATAN(EXP(proj%dlon*(proj%rsw + j-proj%knownj)))*deg_per_rad - 90.
  1083. lon = (i-proj%knowni)*proj%dlon*deg_per_rad + proj%lon1
  1084. IF (lon.GT.180.) lon = lon - 360.
  1085. IF (lon.LT.-180.) lon = lon + 360.
  1086. RETURN
  1087. END SUBROUTINE ijll_merc
  1088. SUBROUTINE llij_latlon(lat, lon, proj, i, j)
  1089. ! Compute the i/j location of a lat/lon on a LATLON grid.
  1090. IMPLICIT NONE
  1091. REAL, INTENT(IN) :: lat
  1092. REAL, INTENT(IN) :: lon
  1093. TYPE(proj_info), INTENT(IN) :: proj
  1094. REAL, INTENT(OUT) :: i
  1095. REAL, INTENT(OUT) :: j
  1096. REAL :: deltalat
  1097. REAL :: deltalon
  1098. ! Compute deltalat and deltalon as the difference between the input
  1099. ! lat/lon and the origin lat/lon
  1100. deltalat = lat - proj%lat1
  1101. deltalon = lon - proj%lon1
  1102. ! Compute i/j
  1103. i = deltalon/proj%loninc
  1104. j = deltalat/proj%latinc
  1105. i = i + proj%knowni
  1106. j = j + proj%knownj
  1107. if ( i < real(proj%nxmin)-0.5 ) i = i + real(proj%nxmax - proj%nxmin + 1)
  1108. if ( i >= real(proj%nxmax)+0.5 ) i = i - real(proj%nxmax - proj%nxmin + 1)
  1109. RETURN
  1110. END SUBROUTINE llij_latlon
  1111. SUBROUTINE ijll_latlon(i, j, proj, lat, lon)
  1112. ! Compute the lat/lon location of an i/j on a LATLON grid.
  1113. IMPLICIT NONE
  1114. REAL, INTENT(IN) :: i
  1115. REAL, INTENT(IN) :: j
  1116. TYPE(proj_info), INTENT(IN) :: proj
  1117. REAL, INTENT(OUT) :: lat
  1118. REAL, INTENT(OUT) :: lon
  1119. REAL :: i_work, j_work
  1120. REAL :: deltalat
  1121. REAL :: deltalon
  1122. i_work = i
  1123. if ( i < real(proj%nxmin)-0.5 ) i_work = i + real(proj%nxmax - proj%nxmin + 1)
  1124. if ( i >= real(proj%nxmax)+0.5 ) i_work = i - real(proj%nxmax - proj%nxmin + 1)
  1125. i_work = i_work - proj%knowni
  1126. j_work = j - proj%knownj
  1127. ! Compute deltalat and deltalon
  1128. deltalat = j_work*proj%latinc
  1129. deltalon = i_work*proj%loninc
  1130. lat = proj%lat1 + deltalat
  1131. lon = proj%lon1 + deltalon
  1132. RETURN
  1133. END SUBROUTINE ijll_latlon
  1134. SUBROUTINE set_cyl(proj)
  1135. implicit none
  1136. ! Arguments
  1137. type(proj_info), intent(inout) :: proj
  1138. proj%hemi = 1.0
  1139. END SUBROUTINE set_cyl
  1140. SUBROUTINE llij_cyl(lat, lon, proj, i, j)
  1141. implicit none
  1142. ! Arguments
  1143. real, intent(in) :: lat, lon
  1144. real, intent(out) :: i, j
  1145. type(proj_info), intent(in) :: proj
  1146. ! Local variables
  1147. real :: deltalat
  1148. real :: deltalon
  1149. ! Compute deltalat and deltalon as the difference between the input
  1150. ! lat/lon and the origin lat/lon
  1151. deltalat = lat - proj%lat1
  1152. ! deltalon = lon - proj%stdlon
  1153. deltalon = lon - proj%lon1
  1154. if (deltalon < 0.) deltalon = deltalon + 360.
  1155. if (deltalon > 360.) deltalon = deltalon - 360.
  1156. ! Compute i/j
  1157. i = deltalon/proj%loninc
  1158. j = deltalat/proj%latinc
  1159. i = i + proj%knowni
  1160. j = j + proj%knownj
  1161. if (i <= 0.) i = i + 360./proj%loninc
  1162. if (i > 360./proj%loninc) i = i - 360./proj%loninc
  1163. END SUBROUTINE llij_cyl
  1164. SUBROUTINE ijll_cyl(i, j, proj, lat, lon)
  1165. implicit none
  1166. ! Arguments
  1167. real, intent(in) :: i, j
  1168. real, intent(out) :: lat, lon
  1169. type(proj_info), intent(in) :: proj
  1170. ! Local variables
  1171. real :: deltalat
  1172. real :: deltalon
  1173. real :: i_work, j_work
  1174. i_work = i - proj%knowni
  1175. j_work = j - proj%knownj
  1176. if (i_work < 0.) i_work = i_work + 360./proj%loninc
  1177. if (i_work >= 360./proj%loninc) i_work = i_work - 360./proj%loninc
  1178. ! Compute deltalat and deltalon
  1179. deltalat = j_work*proj%latinc
  1180. deltalon = i_work*proj%loninc
  1181. lat = deltalat + proj%lat1
  1182. ! lon = deltalon + proj%stdlon
  1183. lon = deltalon + proj%lon1
  1184. if (lon < -180.) lon = lon + 360.
  1185. if (lon > 180.) lon = lon - 360.
  1186. END SUBROUTINE ijll_cyl
  1187. SUBROUTINE set_cassini(proj)
  1188. implicit none
  1189. ! Arguments
  1190. type(proj_info), intent(inout) :: proj
  1191. ! Local variables
  1192. real :: comp_lat, comp_lon
  1193. logical :: global_domain
  1194. proj%hemi = 1.0
  1195. ! Try to determine whether this domain has global coverage
  1196. if (abs(proj%lat1 - proj%latinc/2. + 90.) < 0.001 .and. &
  1197. abs(mod(proj%lon1 - proj%loninc/2. - proj%stdlon,360.)) < 0.001) then
  1198. global_domain = .true.
  1199. else
  1200. global_domain = .false.
  1201. end if
  1202. if (abs(proj%lat0) /= 90. .and. .not.global_domain) then
  1203. call rotate_coords(proj%lat1,proj%lon1,comp_lat,comp_lon,proj%lat0,proj%lon0,proj%stdlon,-1)
  1204. comp_lon = comp_lon + proj%stdlon
  1205. proj%lat1 = comp_lat
  1206. proj%lon1 = comp_lon
  1207. end if
  1208. END SUBROUTINE set_cassini
  1209. SUBROUTINE llij_cassini(lat, lon, proj, i, j)
  1210. implicit none
  1211. ! Arguments
  1212. real, intent(in) :: lat, lon
  1213. real, intent(out) :: i, j
  1214. type(proj_info), intent(in) :: proj
  1215. ! Local variables
  1216. real :: comp_lat, comp_lon
  1217. ! Convert geographic to computational lat/lon
  1218. if ( (abs(proj%lat0) /= 90.) .and. (.not. proj%comp_ll) ) then
  1219. call rotate_coords(lat,lon,comp_lat,comp_lon,proj%lat0,proj%lon0,proj%stdlon,-1)
  1220. comp_lon = comp_lon + proj%stdlon
  1221. else
  1222. comp_lat = lat
  1223. comp_lon = lon
  1224. end if
  1225. ! Convert computational lat/lon to i/j
  1226. call llij_cyl(comp_lat, comp_lon, proj, i, j)
  1227. END SUBROUTINE llij_cassini
  1228. SUBROUTINE ijll_cassini(i, j, proj, lat, lon)
  1229. implicit none
  1230. ! Arguments
  1231. real, intent(in) :: i, j
  1232. real, intent(out) :: lat, lon
  1233. type(proj_info), intent(in) :: proj
  1234. ! Local variables
  1235. real :: comp_lat, comp_lon
  1236. ! Convert i/j to computational lat/lon
  1237. call ijll_cyl(i, j, proj, comp_lat, comp_lon)
  1238. ! Convert computational to geographic lat/lon
  1239. if ( (abs(proj%lat0) /= 90.) .and. (.not. proj%comp_ll) ) then
  1240. comp_lon = comp_lon - proj%stdlon
  1241. call rotate_coords(comp_lat,comp_lon,lat,lon,proj%lat0,proj%lon0,proj%stdlon,1)
  1242. else
  1243. lat = comp_lat
  1244. lon = comp_lon
  1245. end if
  1246. END SUBROUTINE ijll_cassini
  1247. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  1248. ! Purpose: Converts between computational and geographic lat/lon for Cassini
  1249. !
  1250. ! Notes: This routine was provided by Bill Skamarock, 2007-03-27
  1251. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  1252. SUBROUTINE rotate_coords(ilat,ilon,olat,olon,lat_np,lon_np,lon_0,direction)
  1253. IMPLICIT NONE
  1254. REAL, INTENT(IN ) :: ilat, ilon
  1255. REAL, INTENT( OUT) :: olat, olon
  1256. REAL, INTENT(IN ) :: lat_np, lon_np, lon_0
  1257. INTEGER, INTENT(IN ), OPTIONAL :: direction
  1258. ! >=0, default : computational -> geographical
  1259. ! < 0 : geographical -> computational
  1260. REAL :: rlat, rlon
  1261. REAL :: phi_np, lam_np, lam_0, dlam
  1262. REAL :: sinphi, cosphi, coslam, sinlam
  1263. ! Convert all angles to radians
  1264. phi_np = lat_np * rad_per_deg
  1265. lam_np = lon_np * rad_per_deg
  1266. lam_0 = lon_0 * rad_per_deg
  1267. rlat = ilat * rad_per_deg
  1268. rlon = ilon * rad_per_deg
  1269. IF (PRESENT(direction) .AND. (direction < 0)) THEN
  1270. ! The equations are exactly the same except for one small difference
  1271. ! with respect to longitude ...
  1272. dlam = PI - lam_0
  1273. ELSE
  1274. dlam = lam_np
  1275. END IF
  1276. sinphi = COS(phi_np)*COS(rlat)*COS(rlon-dlam) + SIN(phi_np)*SIN(rlat)
  1277. cosphi = SQRT(1.-sinphi*sinphi)
  1278. coslam = SIN(phi_np)*COS(rlat)*COS(rlon-dlam) - COS(phi_np)*SIN(rlat)
  1279. sinlam = COS(rlat)*SIN(rlon-dlam)
  1280. IF ( cosphi /= 0. ) THEN
  1281. coslam = coslam/cosphi
  1282. sinlam = sinlam/cosphi
  1283. END IF
  1284. olat = deg_per_rad*ASIN(sinphi)
  1285. olon = deg_per_rad*(ATAN2(sinlam,coslam)-dlam-lam_0+lam_np)
  1286. ! Both of my F90 text books prefer the DO-EXIT form, and claim it is faster
  1287. ! when optimization is turned on (as we will always do...)
  1288. DO
  1289. IF (olon >= -180.) EXIT
  1290. olon = olon + 360.
  1291. END DO
  1292. DO
  1293. IF (olon <= 180.) EXIT
  1294. olon = olon - 360.
  1295. END DO
  1296. END SUBROUTINE rotate_coords
  1297. SUBROUTINE llij_rotlatlon(lat, lon, proj, i_real, j_real)
  1298. IMPLICIT NONE
  1299. ! Arguments
  1300. REAL, INTENT(IN) :: lat, lon
  1301. REAL :: i, j
  1302. REAL, INTENT(OUT) :: i_real, j_real
  1303. TYPE (proj_info), INTENT(IN) :: proj
  1304. ! Local variables
  1305. INTEGER :: ii,imt,jj,jmt,k,krows,ncol,nrow,iri
  1306. REAL(KIND=HIGH) :: dphd,dlmd !Grid increments, degrees
  1307. REAL(KIND=HIGH) :: glatd !Geographic latitude, positive north
  1308. REAL(KIND=HIGH) :: glond !Geographic longitude, positive west
  1309. REAL(KIND=HIGH) :: col,d1,d2,d2r,dlm,dlm1,dlm2,dph,glat,glon, &
  1310. pi,r2d,row,tlat,tlat1,tlat2, &
  1311. tlon,tlon1,tlon2,tph0,tlm0,x,y,z
  1312. glatd = lat
  1313. glond = -lon
  1314. dphd = proj%phi/REAL((proj%jydim-1)/2)
  1315. dlmd = proj%lambda/REAL(proj%ixdim-1)
  1316. pi = ACOS(-1.0)
  1317. d2r = pi/180.
  1318. r2d = 1./d2r
  1319. imt = 2*proj%ixdim-1
  1320. jmt = proj%jydim/2+1
  1321. glat = glatd*d2r
  1322. glon = glond*d2r
  1323. dph = dphd*d2r
  1324. dlm = dlmd*d2r
  1325. tph0 = proj%lat1*d2r
  1326. tlm0 = -proj%lon1*d2r
  1327. x = COS(tph0)*COS(glat)*COS(glon-tlm0)+SIN(tph0)*SIN(glat)
  1328. y = -COS(glat)*SIN(glon-tlm0)
  1329. z = COS(tph0)*SIN(glat)-SIN(tph0)*COS(glat)*COS(glon-tlm0)
  1330. tlat = r2d*ATAN(z/SQRT(x*x+y*y))
  1331. tlon = r2d*ATAN(y/x)
  1332. row = tlat/dphd+jmt
  1333. col = tlon/dlmd+proj%ixdim
  1334. if ( (row - INT(row)) .gt. 0.999) then
  1335. row = row + 0.0002
  1336. else if ( (col - INT(col)) .gt. 0.999) then
  1337. col = col + 0.0002
  1338. end if
  1339. nrow = INT(row)
  1340. ncol = INT(col)
  1341. ! nrow = NINT(row)
  1342. ! ncol = NINT(col)
  1343. tlat = tlat*d2r
  1344. tlon = tlon*d2r
  1345. IF (proj%stagger == HH) THEN
  1346. IF (mod(nrow,2) .eq. 0) then
  1347. i_real = col / 2.0
  1348. ELSE
  1349. i_real = col / 2.0 + 0.5
  1350. ENDIF
  1351. j_real=row
  1352. IF ((abs(MOD(nrow,2)) == 1 .AND. abs(MOD(ncol,2)) == 1) .OR. &
  1353. (MOD(nrow,2) == 0 .AND. MOD(ncol,2) == 0)) THEN
  1354. tlat1 = (nrow-jmt)*dph
  1355. tlat2 = tlat1+dph
  1356. tlon1 = (ncol-proj%ixdim)*dlm
  1357. tlon2 = tlon1+dlm
  1358. dlm1 = tlon-tlon1
  1359. dlm2 = tlon-tlon2
  1360. d1 = ACOS(COS(tlat)*COS(tlat1)*COS(dlm1)+SIN(tlat)*SIN(tlat1))
  1361. d2 = ACOS(COS(tlat)*COS(tlat2)*COS(dlm2)+SIN(tlat)*SIN(tlat2))
  1362. IF (d1 > d2) THEN
  1363. nrow = nrow+1
  1364. ncol = ncol+1
  1365. END IF
  1366. ELSE
  1367. tlat1 = (nrow+1-jmt)*dph
  1368. tlat2 = tlat1-dph
  1369. tlon1 = (ncol-proj%ixdim)*dlm
  1370. tlon2 = tlon1+dlm
  1371. dlm1 = tlon-tlon1
  1372. dlm2 = tlon-tlon2
  1373. d1 = ACOS(COS(tlat)*COS(tlat1)*COS(dlm1)+SIN(tlat)*SIN(tlat1))
  1374. d2 = ACOS(COS(tlat)*COS(tlat2)*COS(dlm2)+SIN(tlat)*SIN(tlat2))
  1375. IF (d1 < d2) THEN
  1376. nrow = nrow+1
  1377. ELSE
  1378. ncol = ncol+1
  1379. END IF
  1380. END IF
  1381. ELSE IF (proj%stagger == VV) THEN
  1382. IF (mod(nrow,2) .eq. 0) then
  1383. i_real = col / 2.0 + 0.5
  1384. ELSE
  1385. i_real = col / 2.0
  1386. ENDIF
  1387. j_real=row
  1388. IF ((MOD(nrow,2) == 0 .AND. abs(MOD(ncol,2)) == 1) .OR. &
  1389. (abs(MOD(nrow,2)) == 1 .AND. MOD(ncol,2) == 0)) THEN
  1390. tlat1 = (nrow-jmt)*dph
  1391. tlat2 = tlat1+dph
  1392. tlon1 = (ncol-proj%ixdim)*dlm
  1393. tlon2 = tlon1+dlm
  1394. dlm1 = tlon-tlon1
  1395. dlm2 = tlon-tlon2
  1396. d1 = ACOS(COS(tlat)*COS(tlat1)*COS(dlm1)+SIN(tlat)*SIN(tlat1))
  1397. d2 = ACOS(COS(tlat)*COS(tlat2)*COS(dlm2)+SIN(tlat)*SIN(tlat2))
  1398. IF (d1 > d2) THEN
  1399. nrow = nrow+1
  1400. ncol = ncol+1
  1401. END IF
  1402. ELSE
  1403. tlat1 = (nrow+1-jmt)*dph
  1404. tlat2 = tlat1-dph
  1405. tlon1 = (ncol-proj%ixdim)*dlm
  1406. tlon2 = tlon1+dlm
  1407. dlm1 = tlon-tlon1
  1408. dlm2 = tlon-tlon2
  1409. d1 = ACOS(COS(tlat)*COS(tlat1)*COS(dlm1)+SIN(tlat)*SIN(tlat1))
  1410. d2 = ACOS(COS(tlat)*COS(tlat2)*COS(dlm2)+SIN(tlat)*SIN(tlat2))
  1411. IF (d1 < d2) THEN
  1412. nrow = nrow+1
  1413. ELSE
  1414. ncol = ncol+1
  1415. END IF
  1416. END IF
  1417. END IF
  1418. !!! Added next line as a Kludge - not yet understood why needed
  1419. if (ncol .le. 0) ncol=ncol-1
  1420. jj = nrow
  1421. ii = ncol/2
  1422. IF (proj%stagger == HH) THEN
  1423. IF (abs(MOD(jj,2)) == 1) ii = ii+1
  1424. ELSE IF (proj%stagger == VV) THEN
  1425. IF (MOD(jj,2) == 0) ii=ii+1
  1426. END IF
  1427. i = REAL(ii)
  1428. j = REAL(jj)
  1429. END SUBROUTINE llij_rotlatlon
  1430. SUBROUTINE ijll_rotlatlon(i, j, proj, lat,lon)
  1431. IMPLICIT NONE
  1432. ! Arguments
  1433. REAL, INTENT(IN) :: i, j
  1434. REAL, INTENT(OUT) :: lat, lon
  1435. TYPE (proj_info), INTENT(IN) :: proj
  1436. ! Local variables
  1437. INTEGER :: ih,jh
  1438. REAL :: jj
  1439. INTEGER :: midcol,midrow,ncol,iadd1,iadd2,imt,jh2,knrow,krem,kv,nrow
  1440. REAL :: dphd,dlmd !Grid increments, degrees
  1441. REAL(KIND=HIGH) :: arg1,arg2,d2r,fctr,glatr,glatd,glond,pi, &
  1442. r2d,tlatd,tlond,tlatr,tlonr,tlm0,tph0
  1443. REAL :: col
  1444. jj = j
  1445. if ( (j - INT(j)) .gt. 0.999) then
  1446. jj = j + 0.0002
  1447. endif
  1448. jh = INT(jj)
  1449. dphd = proj%phi/REAL((proj%jydim-1)/2)
  1450. dlmd = proj%lambda/REAL(proj%ixdim-1)
  1451. pi = ACOS(-1.0)
  1452. d2r = pi/180.
  1453. r2d = 1./d2r
  1454. tph0 = proj%lat1*d2r
  1455. tlm0 = -proj%lon1*d2r
  1456. midrow = (proj%jydim+1)/2
  1457. midcol = proj%ixdim
  1458. col = 2*i-1+abs(MOD(jh+1,2))
  1459. tlatd = (jj-midrow)*dphd
  1460. tlond = (col-midcol)*dlmd
  1461. IF (proj%stagger == VV) THEN
  1462. if (mod(jh,2) .eq. 0) then
  1463. tlond = tlond - DLMD
  1464. else
  1465. tlond = tlond + DLMD
  1466. end if
  1467. END IF
  1468. tlatr = tlatd*d2r
  1469. tlonr = tlond*d2r
  1470. arg1 = SIN(tlatr)*COS(tph0)+COS(tlatr)*SIN(tph0)*COS(tlonr)
  1471. glatr = ASIN(arg1)
  1472. glatd = glatr*r2d
  1473. arg2 = COS(tlatr)*COS(tlonr)/(COS(glatr)*COS(tph0))-TAN(glatr)*TAN(tph0)
  1474. IF (ABS(arg2) > 1.) arg2 = ABS(arg2)/arg2
  1475. fctr = 1.
  1476. IF (tlond > 0.) fctr = -1.
  1477. glond = tlm0*r2d+fctr*ACOS(arg2)*r2d
  1478. lat = glatd
  1479. lon = -glond
  1480. IF (lon .GT. +180.) lon = lon - 360.
  1481. IF (lon .LT. -180.) lon = lon + 360.
  1482. END SUBROUTINE ijll_rotlatlon
  1483. SUBROUTINE set_gauss(proj)
  1484. IMPLICIT NONE
  1485. ! Argument
  1486. type (proj_info), intent(inout) :: proj
  1487. ! Initialize the array that will hold the Gaussian latitudes.
  1488. IF ( ASSOCIATED( proj%gauss_lat ) ) THEN
  1489. DEALLOCATE ( proj%gauss_lat )
  1490. END IF
  1491. ! Get the needed space for our array.
  1492. ALLOCATE ( proj%gauss_lat(proj%nlat*2) )
  1493. ! Compute the Gaussian latitudes.
  1494. CALL gausll( proj%nlat*2 , proj%gauss_lat )
  1495. ! Now, these could be upside down from what we want, so let's check.
  1496. ! We take advantage of the equatorial symmetry to remove any sort of
  1497. ! array re-ordering.
  1498. IF ( ABS(proj%gauss_lat(1) - proj%lat1) .GT. 0.01 ) THEN
  1499. proj%gauss_lat = -1. * proj%gauss_lat
  1500. END IF
  1501. ! Just a sanity check.
  1502. IF ( ABS(proj%gauss_lat(1) - proj%lat1) .GT. 0.01 ) THEN
  1503. PRINT '(A)','Oops, something is not right with the Gaussian latitude computation.'
  1504. PRINT '(A,F8.3,A)','The input data gave the starting latitude as ',proj%lat1,'.'
  1505. PRINT '(A,F8.3,A)','This routine computed the starting latitude as +-',ABS(proj%gauss_lat(1)),'.'
  1506. PRINT '(A,F8.3,A)','The difference is larger than 0.01 degrees, which is not expected.'
  1507. call mprintf(.true.,ERROR,'Gaussian_latitude_computation')
  1508. END IF
  1509. END SUBROUTINE set_gauss
  1510. SUBROUTINE gausll ( nlat , lat_sp )
  1511. IMPLICIT NONE
  1512. INTEGER :: nlat , i
  1513. REAL (KIND=HIGH) , PARAMETER :: pi = 3.141592653589793
  1514. REAL (KIND=HIGH) , DIMENSION(nlat) :: cosc , gwt , sinc , colat , wos2 , lat
  1515. REAL , DIMENSION(nlat) :: lat_sp
  1516. CALL lggaus(nlat, cosc, gwt, sinc, colat, wos2)
  1517. DO i = 1, nlat
  1518. lat(i) = ACOS(sinc(i)) * 180._HIGH / pi
  1519. IF (i.gt.nlat/2) lat(i) = -lat(i)
  1520. END DO
  1521. lat_sp = REAL(lat)
  1522. END SUBROUTINE gausll
  1523. SUBROUTINE lggaus( nlat, cosc, gwt, sinc, colat, wos2 )
  1524. IMPLICIT NONE
  1525. ! LGGAUS finds the Gaussian latitudes by finding the roots of the
  1526. ! ordinary Legendre polynomial of degree NLAT using Newton's
  1527. ! iteration method.
  1528. ! On entry:
  1529. integer NLAT ! the number of latitudes (degree of the polynomial)
  1530. ! On exit: for each Gaussian latitude
  1531. ! COSC - cos(colatitude) or sin(latitude)
  1532. ! GWT - the Gaussian weights
  1533. ! SINC - sin(colatitude) or cos(latitude)
  1534. ! COLAT - the colatitudes in radians
  1535. ! WOS2 - Gaussian weight over sin**2(colatitude)
  1536. REAL (KIND=HIGH) , DIMENSION(nlat) :: cosc , gwt , sinc , colat , wos2
  1537. REAL (KIND=HIGH) , PARAMETER :: pi = 3.141592653589793
  1538. ! Convergence criterion for iteration of cos latitude
  1539. REAL , PARAMETER :: xlim = 1.0E-14
  1540. INTEGER :: nzero, i, j
  1541. REAL (KIND=HIGH) :: fi, fi1, a, b, g, gm, gp, gt, delta, c, d
  1542. ! The number of zeros between pole and equator
  1543. nzero = nlat/2
  1544. ! Set first guess for cos(colat)
  1545. DO i=1,nzero
  1546. cosc(i) = SIN( (i-0.5)*pi/nlat + pi*0.5 )
  1547. END DO
  1548. ! Constants for determining the derivative of the polynomial
  1549. fi = nlat
  1550. fi1 = fi+1.0
  1551. a = fi*fi1 / SQRT(4.0*fi1*fi1-1.0)
  1552. b = fi1*fi / SQRT(4.0*fi*fi-1.0)
  1553. ! Loop over latitudes, iterating the search for each root
  1554. DO i=1,nzero
  1555. j=0
  1556. ! Determine the value of the ordinary Legendre polynomial for
  1557. ! the current guess root
  1558. DO
  1559. CALL lgord( g, cosc(i), nlat )
  1560. ! Determine the derivative of the polynomial at this point
  1561. CALL lgord( gm, cosc(i), nlat-1 )
  1562. CALL lgord( gp, cosc(i), nlat+1 )
  1563. gt = (cosc(i)*cosc(i)-1.0) / (a*gp-b*gm)
  1564. ! Update the estimate of the root
  1565. delta = g*gt
  1566. cosc(i) = cosc(i) - delta
  1567. ! If convergence criterion has not been met, keep trying
  1568. j = j+1
  1569. IF( ABS(delta).GT.xlim ) CYCLE
  1570. ! Determine the Gaussian weights
  1571. c = 2.0 *( 1.0-cosc(i)*cosc(i) )
  1572. CALL lgord( d, cosc(i), nlat-1 )
  1573. d = d*d*fi*fi
  1574. gwt(i) = c *( fi-0.5 ) / d
  1575. EXIT
  1576. END DO
  1577. END DO
  1578. ! Determine the colatitudes and sin(colat) and weights over sin**2
  1579. DO i=1,nzero
  1580. colat(i)= ACOS(cosc(i))
  1581. sinc(i) = SIN(colat(i))
  1582. wos2(i) = gwt(i) /( sinc(i)*sinc(i) )
  1583. END DO
  1584. ! If NLAT is odd, set values at the equator
  1585. IF( MOD(nlat,2) .NE. 0 ) THEN
  1586. i = nzero+1
  1587. cosc(i) = 0.0
  1588. c = 2.0
  1589. CALL lgord( d, cosc(i), nlat-1 )
  1590. d = d*d*fi*fi
  1591. gwt(i) = c *( fi-0.5 ) / d
  1592. colat(i)= pi*0.5
  1593. sinc(i) = 1.0
  1594. wos2(i) = gwt(i)
  1595. END IF
  1596. ! Determine the southern hemisphere values by symmetry
  1597. DO i=nlat-nzero+1,nlat
  1598. cosc(i) =-cosc(nlat+1-i)
  1599. gwt(i) = gwt(nlat+1-i)
  1600. colat(i)= pi-colat(nlat+1-i)
  1601. sinc(i) = sinc(nlat+1-i)
  1602. wos2(i) = wos2(nlat+1-i)
  1603. END DO
  1604. END SUBROUTINE lggaus
  1605. SUBROUTINE lgord( f, cosc, n )
  1606. IMPLICIT NONE
  1607. ! LGORD calculates the value of an ordinary Legendre polynomial at a
  1608. ! specific latitude.
  1609. ! On entry:
  1610. ! cosc - COS(colatitude)
  1611. ! n - the degree of the polynomial
  1612. ! On exit:
  1613. ! f - the value of the Legendre polynomial of degree N at
  1614. ! latitude ASIN(cosc)
  1615. REAL (KIND=HIGH) :: s1, c4, a, b, fk, f, cosc, colat, c1, fn, ang
  1616. INTEGER :: n, k
  1617. ! Determine the colatitude
  1618. colat = ACOS(cosc)
  1619. c1 = SQRT(2.0_HIGH)
  1620. DO k=1,n
  1621. c1 = c1 * SQRT( 1.0 - 1.0/(4*k*k) )
  1622. END DO
  1623. fn = n
  1624. ang= fn * colat
  1625. s1 = 0.0
  1626. c4 = 1.0
  1627. a =-1.0
  1628. b = 0.0
  1629. DO k=0,n,2
  1630. IF (k.eq.n) c4 = 0.5 * c4
  1631. s1 = s1 + c4 * COS(ang)
  1632. a = a + 2.0
  1633. b = b + 1.0
  1634. fk = k
  1635. ang= colat * (fn-fk-2.0)
  1636. c4 = ( a * (fn-b+1.0) / ( b * (fn+fn-a) ) ) * c4
  1637. END DO
  1638. f = s1 * c1
  1639. END SUBROUTINE lgord
  1640. SUBROUTINE llij_gauss (lat, lon, proj, i, j)
  1641. IMPLICIT NONE
  1642. REAL , INTENT(IN) :: lat, lon
  1643. REAL , INTENT(OUT) :: i, j
  1644. TYPE (proj_info), INTENT(IN) :: proj
  1645. INTEGER :: n , n_low
  1646. LOGICAL :: found = .FALSE.
  1647. REAL :: diff_1 , diff_nlat
  1648. ! The easy one first, get the x location. The calling routine has already made
  1649. ! sure that the necessary assumptions concerning the sign of the deltalon and the
  1650. ! relative east/west'ness of the longitude and the starting longitude are consistent
  1651. ! to allow this easy computation.
  1652. i = ( lon - proj%lon1 ) / proj%loninc + 1.
  1653. ! Since this is a global data set, we need to be concerned about wrapping the
  1654. ! fields around the globe.
  1655. ! IF ( ( proj%loninc .GT. 0 ) .AND. &
  1656. ! ( FLOOR((lon-proj%lon1)/proj%loninc) + 1 .GE. proj%ixdim ) .AND. &
  1657. ! ( lon + proj%loninc .GE. proj%lon1 + 360 ) ) THEN
  1658. !! BUG: We may need to set proj%wrap, but proj is intent(in)
  1659. !! WHAT IS THIS USED FOR?
  1660. !! proj%wrap = .TRUE.
  1661. ! i = proj%ixdim
  1662. ! ELSE IF ( ( proj%loninc .LT. 0 ) .AND. &
  1663. ! ( FLOOR((lon-proj%lon1)/proj%loninc) + 1 .GE. proj%ixdim ) .AND. &
  1664. ! ( lon + proj%loninc .LE. proj%lon1 - 360 ) ) THEN
  1665. ! ! BUG: We may need to set proj%wrap, but proj is intent(in)
  1666. ! ! WHAT IS THIS USED FOR?
  1667. ! ! proj%wrap = .TRUE.
  1668. ! i = proj%ixdim
  1669. ! END IF
  1670. ! Yet another quicky test, can we find bounding values? If not, then we may be
  1671. ! dealing with putting data to a polar projection, so just give them them maximal
  1672. ! value for the location. This is an OK assumption for the interpolation across the
  1673. ! top of the pole, given how close the longitude lines are.
  1674. IF ( ABS(lat) .GT. ABS(proj%gauss_lat(1)) ) THEN
  1675. diff_1 = lat - proj%gauss_lat(1)
  1676. diff_nlat = lat - proj%gauss_lat(proj%nlat*2)
  1677. IF ( ABS(diff_1) .LT. ABS(diff_nlat) ) THEN
  1678. j = 1
  1679. ELSE
  1680. j = proj%nlat*2
  1681. END IF
  1682. ! If the latitude is between the two bounding values, we have to search and interpolate.
  1683. ELSE
  1684. DO n = 1 , proj%nlat*2 -1
  1685. IF ( ( proj%gauss_lat(n) - lat ) * ( proj%gauss_lat(n+1) - lat ) .LE. 0 ) THEN
  1686. found = .TRUE.
  1687. n_low = n
  1688. EXIT
  1689. END IF
  1690. END DO
  1691. ! Everything still OK?
  1692. IF ( .NOT. found ) THEN
  1693. PRINT '(A)','Troubles in river city. No bounding values of latitude found in the Gaussian routines.'
  1694. call mprintf(.true.,ERROR,'Gee_no_bounding_lats_Gaussian')
  1695. END IF
  1696. j = ( ( proj%gauss_lat(n_low) - lat ) * ( n_low + 1 ) + &
  1697. ( lat - proj%gauss_lat(n_low+1) ) * ( n_low ) ) / &
  1698. ( proj%gauss_lat(n_low) - proj%gauss_lat(n_low+1) )
  1699. END IF
  1700. if ( i < real(proj%nxmin)-0.5 ) i = i + real(proj%nxmax - proj%nxmin + 1)
  1701. if ( i >= real(proj%nxmax)+0.5 ) i = i - real(proj%nxmax - proj%nxmin + 1)
  1702. END SUBROUTINE llij_gauss
  1703. END MODULE map_utils