remap_bicubic.f 28 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844
  1. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  2. !
  3. ! this module contains necessary routines for performing an
  4. ! bicubic interpolation.
  5. !
  6. !-----------------------------------------------------------------------
  7. !
  8. ! CVS:$Id: remap_bicubic.f,v 1.5 2001/08/22 18:20:41 pwjones Exp $
  9. !
  10. ! Copyright (c) 1997, 1998 the Regents of the University of
  11. ! California.
  12. !
  13. ! This software and ancillary information (herein called software)
  14. ! called SCRIP is made available under the terms described here.
  15. ! The software has been approved for release with associated
  16. ! LA-CC Number 98-45.
  17. !
  18. ! Unless otherwise indicated, this software has been authored
  19. ! by an employee or employees of the University of California,
  20. ! operator of the Los Alamos National Laboratory under Contract
  21. ! No. W-7405-ENG-36 with the U.S. Department of Energy. The U.S.
  22. ! Government has rights to use, reproduce, and distribute this
  23. ! software. The public may copy and use this software without
  24. ! charge, provided that this Notice and any statement of authorship
  25. ! are reproduced on all copies. Neither the Government nor the
  26. ! University makes any warranty, express or implied, or assumes
  27. ! any liability or responsibility for the use of this software.
  28. !
  29. ! If software is modified to produce derivative works, such modified
  30. ! software should be clearly marked, so as not to confuse it with
  31. ! the version available from Los Alamos National Laboratory.
  32. !
  33. !***********************************************************************
  34. module remap_bicubic
  35. !-----------------------------------------------------------------------
  36. use kinds_mod ! defines common data types
  37. use constants ! defines common constants
  38. use grids ! module containing grid info
  39. use remap_vars ! module containing remap info
  40. implicit none
  41. !-----------------------------------------------------------------------
  42. integer (kind=int_kind), parameter ::
  43. & max_iter = 100 ! max iteration count for i,j iteration
  44. real (kind=dbl_kind), parameter ::
  45. & converge = 1.e-10_dbl_kind ! convergence criterion
  46. !***********************************************************************
  47. contains
  48. !***********************************************************************
  49. subroutine remap_bicub
  50. !-----------------------------------------------------------------------
  51. !
  52. ! this routine computes the weights for a bicubic interpolation.
  53. !
  54. !-----------------------------------------------------------------------
  55. !-----------------------------------------------------------------------
  56. !
  57. ! local variables
  58. !
  59. !-----------------------------------------------------------------------
  60. integer (kind=int_kind) :: n,icount,
  61. & dst_add, ! destination address
  62. & iter, ! iteration counter
  63. & nmap ! index of current map being computed
  64. integer (kind=int_kind), dimension(4) ::
  65. & src_add ! address for the four source points
  66. real (kind=dbl_kind), dimension(4) ::
  67. & src_lats, ! latitudes of four bilinear corners
  68. & src_lons ! longitudes of four bilinear corners
  69. real (kind=dbl_kind), dimension(4,4) ::
  70. & wgts ! bicubic weights for four corners
  71. real (kind=dbl_kind) ::
  72. & plat, plon, ! lat/lon coords of destination point
  73. & iguess, jguess, ! current guess for bilinear coordinate
  74. & thguess, phguess, ! current guess for lat/lon coordinate
  75. & deli, delj, ! corrections to i,j
  76. & dth1, dth2, dth3, ! some latitude differences
  77. & dph1, dph2, dph3, ! some longitude differences
  78. & dthp, dphp, ! difference between point and sw corner
  79. & mat1, mat2, mat3, mat4, ! matrix elements
  80. & determinant, ! matrix determinant
  81. & sum_wgts, ! sum of weights for normalization
  82. & w1,w2,w3,w4,w5,w6,w7,w8, ! 16 bicubic weight functions
  83. & w9,w10,w11,w12,w13,w14,w15,w16
  84. !-----------------------------------------------------------------------
  85. !
  86. ! compute mappings from grid1 to grid2
  87. !
  88. !-----------------------------------------------------------------------
  89. nmap = 1
  90. if (grid1_rank /= 2) then
  91. stop 'Can not do bicubic interpolation when grid_rank /= 2'
  92. endif
  93. !***
  94. !*** loop over destination grid
  95. !***
  96. grid_loop1: do dst_add = 1, grid2_size
  97. if (.not. grid2_mask(dst_add)) cycle grid_loop1
  98. plat = grid2_center_lat(dst_add)
  99. plon = grid2_center_lon(dst_add)
  100. !-----------------------------------------------------------------------
  101. !
  102. ! find nearest square of grid points on source grid
  103. !
  104. !-----------------------------------------------------------------------
  105. call grid_search_bicub(src_add, src_lats, src_lons,
  106. & plat, plon, grid1_dims,
  107. & grid1_center_lat, grid1_center_lon,
  108. & grid1_bound_box, bin_addr1, bin_addr2)
  109. !***
  110. !*** check to see if points are land points
  111. !***
  112. if (src_add(1) > 0) then
  113. do n=1,4
  114. if (.not. grid1_mask(src_add(n))) src_add(1) = 0
  115. end do
  116. endif
  117. !-----------------------------------------------------------------------
  118. !
  119. ! if point found, find local i,j coordinates for weights
  120. !
  121. !-----------------------------------------------------------------------
  122. if (src_add(1) > 0) then
  123. grid2_frac(dst_add) = one
  124. !***
  125. !*** iterate to find i,j for bicubic approximation
  126. !***
  127. dth1 = src_lats(2) - src_lats(1)
  128. dth2 = src_lats(4) - src_lats(1)
  129. dth3 = src_lats(3) - src_lats(2) - dth2
  130. dph1 = src_lons(2) - src_lons(1)
  131. dph2 = src_lons(4) - src_lons(1)
  132. dph3 = src_lons(3) - src_lons(2)
  133. if (dph1 > three*pih) dph1 = dph1 - pi2
  134. if (dph2 > three*pih) dph2 = dph2 - pi2
  135. if (dph3 > three*pih) dph3 = dph3 - pi2
  136. if (dph1 < -three*pih) dph1 = dph1 + pi2
  137. if (dph2 < -three*pih) dph2 = dph2 + pi2
  138. if (dph3 < -three*pih) dph3 = dph3 + pi2
  139. dph3 = dph3 - dph2
  140. iguess = half
  141. jguess = half
  142. iter_loop1: do iter=1,max_iter
  143. dthp = plat - src_lats(1) - dth1*iguess -
  144. & dth2*jguess - dth3*iguess*jguess
  145. dphp = plon - src_lons(1)
  146. if (dphp > three*pih) dphp = dphp - pi2
  147. if (dphp < -three*pih) dphp = dphp + pi2
  148. dphp = dphp - dph1*iguess - dph2*jguess -
  149. & dph3*iguess*jguess
  150. mat1 = dth1 + dth3*jguess
  151. mat2 = dth2 + dth3*iguess
  152. mat3 = dph1 + dph3*jguess
  153. mat4 = dph2 + dph3*iguess
  154. determinant = mat1*mat4 - mat2*mat3
  155. deli = (dthp*mat4 - mat2*dphp)/determinant
  156. delj = (mat1*dphp - dthp*mat3)/determinant
  157. if (abs(deli) < converge .and.
  158. & abs(delj) < converge) exit iter_loop1
  159. iguess = iguess + deli
  160. jguess = jguess + delj
  161. end do iter_loop1
  162. if (iter <= max_iter) then
  163. !-----------------------------------------------------------------------
  164. !
  165. ! successfully found i,j - compute weights
  166. !
  167. !-----------------------------------------------------------------------
  168. wgts(1,1) = (one - jguess**2*(three-two*jguess))*
  169. & (one - iguess**2*(three-two*iguess))
  170. wgts(1,2) = (one - jguess**2*(three-two*jguess))*
  171. & iguess**2*(three-two*iguess)
  172. wgts(1,3) = jguess**2*(three-two*jguess)*
  173. & iguess**2*(three-two*iguess)
  174. wgts(1,4) = jguess**2*(three-two*jguess)*
  175. & (one - iguess**2*(three-two*iguess))
  176. wgts(2,1) = (one - jguess**2*(three-two*jguess))*
  177. & iguess*(iguess-one)**2
  178. wgts(2,2) = (one - jguess**2*(three-two*jguess))*
  179. & iguess**2*(iguess-one)
  180. wgts(2,3) = jguess**2*(three-two*jguess)*
  181. & iguess**2*(iguess-one)
  182. wgts(2,4) = jguess**2*(three-two*jguess)*
  183. & iguess*(iguess-one)**2
  184. wgts(3,1) = jguess*(jguess-one)**2*
  185. & (one - iguess**2*(three-two*iguess))
  186. wgts(3,2) = jguess*(jguess-one)**2*
  187. & iguess**2*(three-two*iguess)
  188. wgts(3,3) = jguess**2*(jguess-one)*
  189. & iguess**2*(three-two*iguess)
  190. wgts(3,4) = jguess**2*(jguess-one)*
  191. & (one - iguess**2*(three-two*iguess))
  192. wgts(4,1) = iguess*(iguess-one)**2*
  193. & jguess*(jguess-one)**2
  194. wgts(4,2) = iguess**2*(iguess-one)*
  195. & jguess*(jguess-one)**2
  196. wgts(4,3) = iguess**2*(iguess-one)*
  197. & jguess**2*(jguess-one)
  198. wgts(4,4) = iguess*(iguess-one)**2*
  199. & jguess**2*(jguess-one)
  200. call store_link_bicub(dst_add, src_add, wgts, nmap)
  201. else
  202. stop 'Iteration for i,j exceed max iteration count'
  203. endif
  204. !-----------------------------------------------------------------------
  205. !
  206. ! search for bilinear failed - use a distance-weighted
  207. ! average instead (this is typically near the pole)
  208. !
  209. !-----------------------------------------------------------------------
  210. else if (src_add(1) < 0) then
  211. src_add = abs(src_add)
  212. icount = 0
  213. do n=1,4
  214. if (grid1_mask(src_add(n))) then
  215. icount = icount + 1
  216. else
  217. src_lats(n) = zero
  218. endif
  219. end do
  220. if (icount > 0) then
  221. !*** renormalize weights
  222. sum_wgts = sum(src_lats)
  223. wgts(1,1) = src_lats(1)/sum_wgts
  224. wgts(1,2) = src_lats(2)/sum_wgts
  225. wgts(1,3) = src_lats(3)/sum_wgts
  226. wgts(1,4) = src_lats(4)/sum_wgts
  227. wgts(2:4,:) = zero
  228. grid2_frac(dst_add) = one
  229. call store_link_bicub(dst_add, src_add, wgts, nmap)
  230. endif
  231. endif
  232. end do grid_loop1
  233. !-----------------------------------------------------------------------
  234. !
  235. ! compute mappings from grid2 to grid1 if necessary
  236. !
  237. !-----------------------------------------------------------------------
  238. if (num_maps > 1) then
  239. nmap = 2
  240. if (grid2_rank /= 2) then
  241. stop 'Can not do bicubic interpolation when grid_rank /= 2'
  242. endif
  243. !***
  244. !*** loop over destination grid
  245. !***
  246. grid_loop2: do dst_add = 1, grid1_size
  247. if (.not. grid1_mask(dst_add)) cycle grid_loop2
  248. plat = grid1_center_lat(dst_add)
  249. plon = grid1_center_lon(dst_add)
  250. !***
  251. !*** find nearest square of grid points on source grid
  252. !***
  253. call grid_search_bicub(src_add, src_lats, src_lons,
  254. & plat, plon, grid2_dims,
  255. & grid2_center_lat, grid2_center_lon,
  256. & grid2_bound_box, bin_addr2, bin_addr1)
  257. !***
  258. !*** check to see if points are land points
  259. !***
  260. if (src_add(1) > 0) then
  261. do n=1,4
  262. if (.not. grid2_mask(src_add(n))) src_add(1) = 0
  263. end do
  264. endif
  265. !***
  266. !*** if point found, find i,j coordinates for weights
  267. !***
  268. if (src_add(1) > 0) then
  269. grid1_frac(dst_add) = one
  270. !***
  271. !*** iterate to find i,j for bilinear approximation
  272. !***
  273. dth1 = src_lats(2) - src_lats(1)
  274. dth2 = src_lats(4) - src_lats(1)
  275. dth3 = src_lats(3) - src_lats(2) - dth2
  276. dph1 = src_lons(2) - src_lons(1)
  277. dph2 = src_lons(4) - src_lons(1)
  278. dph3 = src_lons(3) - src_lons(2)
  279. if (dph1 > pi) dph1 = dph1 - pi2
  280. if (dph2 > pi) dph2 = dph2 - pi2
  281. if (dph3 > pi) dph3 = dph3 - pi2
  282. if (dph1 < -pi) dph1 = dph1 + pi2
  283. if (dph2 < -pi) dph2 = dph2 + pi2
  284. if (dph3 < -pi) dph3 = dph3 + pi2
  285. dph3 = dph3 - dph2
  286. iguess = zero
  287. jguess = zero
  288. iter_loop2: do iter=1,max_iter
  289. dthp = plat - src_lats(1) - dth1*iguess -
  290. & dth2*jguess - dth3*iguess*jguess
  291. dphp = plon - src_lons(1)
  292. if (dphp > pi) dphp = dphp - pi2
  293. if (dphp < -pi) dphp = dphp + pi2
  294. dphp = dphp - dph1*iguess - dph2*jguess -
  295. & dph3*iguess*jguess
  296. mat1 = dth1 + dth3*jguess
  297. mat2 = dth2 + dth3*iguess
  298. mat3 = dph1 + dph3*jguess
  299. mat4 = dph2 + dph3*iguess
  300. determinant = mat1*mat4 - mat2*mat3
  301. deli = (dthp*mat4 - mat2*dphp)/determinant
  302. delj = (mat1*dphp - dthp*mat3)/determinant
  303. if (abs(deli) < converge .and.
  304. & abs(delj) < converge) exit iter_loop2
  305. iguess = iguess + deli
  306. jguess = jguess + delj
  307. end do iter_loop2
  308. if (iter <= max_iter) then
  309. !***
  310. !*** successfully found i,j - compute weights
  311. !***
  312. wgts(1,1) = (one - jguess**2*(three-two*jguess))*
  313. & (one - iguess**2*(three-two*iguess))
  314. wgts(1,2) = (one - jguess**2*(three-two*jguess))*
  315. & iguess**2*(three-two*iguess)
  316. wgts(1,3) = jguess**2*(three-two*jguess)*
  317. & iguess**2*(three-two*iguess)
  318. wgts(1,4) = jguess**2*(three-two*jguess)*
  319. & (one - iguess**2*(three-two*iguess))
  320. wgts(2,1) = (one - jguess**2*(three-two*jguess))*
  321. & iguess*(iguess-one)**2
  322. wgts(2,2) = (one - jguess**2*(three-two*jguess))*
  323. & iguess**2*(iguess-one)
  324. wgts(2,3) = jguess**2*(three-two*jguess)*
  325. & iguess**2*(iguess-one)
  326. wgts(2,4) = jguess**2*(three-two*jguess)*
  327. & iguess*(iguess-one)**2
  328. wgts(3,1) = jguess*(jguess-one)**2*
  329. & (one - iguess**2*(three-two*iguess))
  330. wgts(3,2) = jguess*(jguess-one)**2*
  331. & iguess**2*(three-two*iguess)
  332. wgts(3,3) = jguess**2*(jguess-one)*
  333. & iguess**2*(three-two*iguess)
  334. wgts(3,4) = jguess**2*(jguess-one)*
  335. & (one - iguess**2*(three-two*iguess))
  336. wgts(4,1) = iguess*(iguess-one)**2*
  337. & jguess*(jguess-one)**2
  338. wgts(4,2) = iguess**2*(iguess-one)*
  339. & jguess*(jguess-one)**2
  340. wgts(4,3) = iguess**2*(iguess-one)*
  341. & jguess**2*(jguess-one)
  342. wgts(4,4) = iguess*(iguess-one)**2*
  343. & jguess**2*(jguess-one)
  344. call store_link_bicub(dst_add, src_add, wgts, nmap)
  345. else
  346. stop 'Iteration for i,j exceed max iteration count'
  347. endif
  348. !***
  349. !*** search for bilinear failed - us a distance-weighted
  350. !*** average instead
  351. !***
  352. else if (src_add(1) < 0) then
  353. src_add = abs(src_add)
  354. icount = 0
  355. do n=1,4
  356. if (grid2_mask(src_add(n))) then
  357. icount = icount + 1
  358. else
  359. src_lats(n) = zero
  360. endif
  361. end do
  362. if (icount > 0) then
  363. !*** renormalize weights
  364. sum_wgts = sum(src_lats)
  365. wgts(1,1) = src_lats(1)/sum_wgts
  366. wgts(1,2) = src_lats(2)/sum_wgts
  367. wgts(1,3) = src_lats(3)/sum_wgts
  368. wgts(1,4) = src_lats(4)/sum_wgts
  369. wgts(2:4,:) = zero
  370. grid1_frac(dst_add) = one
  371. call store_link_bicub(dst_add, src_add, wgts, nmap)
  372. endif
  373. endif
  374. end do grid_loop2
  375. endif ! nmap=2
  376. !-----------------------------------------------------------------------
  377. end subroutine remap_bicub
  378. !***********************************************************************
  379. subroutine grid_search_bicub(src_add, src_lats, src_lons,
  380. & plat, plon, src_grid_dims,
  381. & src_center_lat, src_center_lon,
  382. & src_bound_box,
  383. & src_bin_add, dst_bin_add)
  384. !-----------------------------------------------------------------------
  385. !
  386. ! this routine finds the location of the search point plat, plon
  387. ! in the source grid and returns the corners needed for a bicubic
  388. ! interpolation.
  389. !
  390. !-----------------------------------------------------------------------
  391. !-----------------------------------------------------------------------
  392. !
  393. ! output variables
  394. !
  395. !-----------------------------------------------------------------------
  396. integer (kind=int_kind), dimension(4), intent(out) ::
  397. & src_add ! address of each corner point enclosing P
  398. real (kind=dbl_kind), dimension(4), intent(out) ::
  399. & src_lats, ! latitudes of the four corner points
  400. & src_lons ! longitudes of the four corner points
  401. !-----------------------------------------------------------------------
  402. !
  403. ! input variables
  404. !
  405. !-----------------------------------------------------------------------
  406. real (kind=dbl_kind), intent(in) ::
  407. & plat, ! latitude of the search point
  408. & plon ! longitude of the search point
  409. integer (kind=int_kind), dimension(2), intent(in) ::
  410. & src_grid_dims ! size of each src grid dimension
  411. real (kind=dbl_kind), dimension(:), intent(in) ::
  412. & src_center_lat, ! latitude of each src grid center
  413. & src_center_lon ! longitude of each src grid center
  414. real (kind=dbl_kind), dimension(:,:), intent(in) ::
  415. & src_bound_box ! bounding box for src grid search
  416. integer (kind=int_kind), dimension(:,:), intent(in) ::
  417. & src_bin_add, ! search bins for restricting
  418. & dst_bin_add ! searches
  419. !-----------------------------------------------------------------------
  420. !
  421. ! local variables
  422. !
  423. !-----------------------------------------------------------------------
  424. integer (kind=int_kind) :: n, next_n, srch_add, ! dummy indices
  425. & nx, ny, ! dimensions of src grid
  426. & min_add, max_add, ! addresses for restricting search
  427. & i, j, jp1, ip1, n_add, e_add, ne_add ! addresses
  428. real (kind=dbl_kind) :: ! vectors for cross-product check
  429. & vec1_lat, vec1_lon,
  430. & vec2_lat, vec2_lon, cross_product, cross_product_last,
  431. & coslat_dst, sinlat_dst, coslon_dst, sinlon_dst,
  432. & dist_min, distance ! for computing dist-weighted avg
  433. !-----------------------------------------------------------------------
  434. !
  435. ! restrict search first using search bins.
  436. !
  437. !-----------------------------------------------------------------------
  438. src_add = 0
  439. min_add = size(src_center_lat)
  440. max_add = 1
  441. do n=1,num_srch_bins
  442. if (plat >= bin_lats(1,n) .and. plat <= bin_lats(2,n) .and.
  443. & plon >= bin_lons(1,n) .and. plon <= bin_lons(2,n)) then
  444. min_add = min(min_add, src_bin_add(1,n))
  445. max_add = max(max_add, src_bin_add(2,n))
  446. endif
  447. end do
  448. !-----------------------------------------------------------------------
  449. !
  450. ! now perform a more detailed search
  451. !
  452. !-----------------------------------------------------------------------
  453. nx = src_grid_dims(1)
  454. ny = src_grid_dims(2)
  455. srch_loop: do srch_add = min_add,max_add
  456. if (plat <= src_bound_box(2,srch_add) .and.
  457. & plat >= src_bound_box(1,srch_add) .and.
  458. & plon <= src_bound_box(4,srch_add) .and.
  459. & plon >= src_bound_box(3,srch_add)) then
  460. !***
  461. !*** we are within bounding box so get really serious
  462. !***
  463. !*** find N,S and NE points to this grid point
  464. j = (srch_add - 1)/nx +1
  465. i = srch_add - (j-1)*nx
  466. if (i < nx) then
  467. ip1 = i + 1
  468. else
  469. ip1 = 1
  470. endif
  471. if (j < ny) then
  472. jp1 = j+1
  473. else
  474. jp1 = 1
  475. endif
  476. n_add = (jp1 - 1)*nx + i
  477. e_add = (j - 1)*nx + ip1
  478. ne_add = (jp1 - 1)*nx + ip1
  479. !***
  480. !*** find N,S and NE lat/lon coords and check bounding box
  481. !***
  482. src_lats(1) = src_center_lat(srch_add)
  483. src_lats(2) = src_center_lat(e_add)
  484. src_lats(3) = src_center_lat(ne_add)
  485. src_lats(4) = src_center_lat(n_add)
  486. src_lons(1) = src_center_lon(srch_add)
  487. src_lons(2) = src_center_lon(e_add)
  488. src_lons(3) = src_center_lon(ne_add)
  489. src_lons(4) = src_center_lon(n_add)
  490. !***
  491. !*** for consistency, we must make sure all lons are in
  492. !*** same 2pi interval
  493. !***
  494. vec1_lon = src_lons(1) - plon
  495. if (vec1_lon > pi) then
  496. src_lons(1) = src_lons(1) - pi2
  497. else if (vec1_lon < -pi) then
  498. src_lons(1) = src_lons(1) + pi2
  499. endif
  500. do n=2,4
  501. vec1_lon = src_lons(n) - src_lons(1)
  502. if (vec1_lon > pi) then
  503. src_lons(n) = src_lons(n) - pi2
  504. else if (vec1_lon < -pi) then
  505. src_lons(n) = src_lons(n) + pi2
  506. endif
  507. end do
  508. corner_loop: do n=1,4
  509. next_n = MOD(n,4) + 1
  510. !***
  511. !*** here we take the cross product of the vector making
  512. !*** up each box side with the vector formed by the vertex
  513. !*** and search point. if all the cross products are
  514. !*** same sign, the point is contained in the box.
  515. !***
  516. vec1_lat = src_lats(next_n) - src_lats(n)
  517. vec1_lon = src_lons(next_n) - src_lons(n)
  518. vec2_lat = plat - src_lats(n)
  519. vec2_lon = plon - src_lons(n)
  520. !***
  521. !*** check for 0,2pi crossings
  522. !***
  523. if (vec1_lon > three*pih) then
  524. vec1_lon = vec1_lon - pi2
  525. else if (vec1_lon < -three*pih) then
  526. vec1_lon = vec1_lon + pi2
  527. endif
  528. if (vec2_lon > three*pih) then
  529. vec2_lon = vec2_lon - pi2
  530. else if (vec2_lon < -three*pih) then
  531. vec2_lon = vec2_lon + pi2
  532. endif
  533. cross_product = vec1_lon*vec2_lat - vec2_lon*vec1_lat
  534. !***
  535. !*** if cross product is less than zero, this cell
  536. !*** doesn't work
  537. !***
  538. if (n==1) cross_product_last = cross_product
  539. if (cross_product*cross_product_last < zero) then
  540. exit corner_loop
  541. else
  542. cross_product_last = cross_product
  543. endif
  544. end do corner_loop
  545. !***
  546. !*** if cross products all positive, we found the location
  547. !***
  548. if (n > 4) then
  549. src_add(1) = srch_add
  550. src_add(2) = e_add
  551. src_add(3) = ne_add
  552. src_add(4) = n_add
  553. return
  554. endif
  555. !***
  556. !*** otherwise move on to next cell
  557. !***
  558. endif !bounding box check
  559. end do srch_loop
  560. !***
  561. !*** if no cell found, point is likely either in a box that
  562. !*** straddles either pole or is outside the grid. fall back
  563. !*** to a distance-weighted average of the four closest
  564. !*** points. go ahead and compute weights here, but store
  565. !*** in src_lats and return -add to prevent the parent
  566. !*** routine from computing bilinear weights
  567. !***
  568. coslat_dst = cos(plat)
  569. sinlat_dst = sin(plat)
  570. coslon_dst = cos(plon)
  571. sinlon_dst = sin(plon)
  572. dist_min = bignum
  573. src_lats = bignum
  574. do srch_add = min_add,max_add
  575. distance = acos(coslat_dst*cos(src_center_lat(srch_add))*
  576. & (coslon_dst*cos(src_center_lon(srch_add)) +
  577. & sinlon_dst*sin(src_center_lon(srch_add)))+
  578. & sinlat_dst*sin(src_center_lat(srch_add)))
  579. if (distance < dist_min) then
  580. sort_loop: do n=1,4
  581. if (distance < src_lats(n)) then
  582. do i=4,n+1,-1
  583. src_add (i) = src_add (i-1)
  584. src_lats(i) = src_lats(i-1)
  585. end do
  586. src_add (n) = -srch_add
  587. src_lats(n) = distance
  588. dist_min = src_lats(4)
  589. exit sort_loop
  590. endif
  591. end do sort_loop
  592. endif
  593. end do
  594. src_lons = one/(src_lats + tiny)
  595. distance = sum(src_lons)
  596. src_lats = src_lons/distance
  597. !-----------------------------------------------------------------------
  598. end subroutine grid_search_bicub
  599. !***********************************************************************
  600. subroutine store_link_bicub(dst_add, src_add, weights, nmap)
  601. !-----------------------------------------------------------------------
  602. !
  603. ! this routine stores the address and weight for four links
  604. ! associated with one destination point in the appropriate address
  605. ! and weight arrays and resizes those arrays if necessary.
  606. !
  607. !-----------------------------------------------------------------------
  608. !-----------------------------------------------------------------------
  609. !
  610. ! input variables
  611. !
  612. !-----------------------------------------------------------------------
  613. integer (kind=int_kind), intent(in) ::
  614. & dst_add, ! address on destination grid
  615. & nmap ! identifies which direction for mapping
  616. integer (kind=int_kind), dimension(4), intent(in) ::
  617. & src_add ! addresses on source grid
  618. real (kind=dbl_kind), dimension(4,4), intent(in) ::
  619. & weights ! array of remapping weights for these links
  620. !-----------------------------------------------------------------------
  621. !
  622. ! local variables
  623. !
  624. !-----------------------------------------------------------------------
  625. integer (kind=int_kind) :: n, ! dummy index
  626. & num_links_old ! placeholder for old link number
  627. !-----------------------------------------------------------------------
  628. !
  629. ! increment number of links and check to see if remap arrays need
  630. ! to be increased to accomodate the new link. then store the
  631. ! link.
  632. !
  633. !-----------------------------------------------------------------------
  634. select case (nmap)
  635. case(1)
  636. num_links_old = num_links_map1
  637. num_links_map1 = num_links_old + 4
  638. if (num_links_map1 > max_links_map1)
  639. & call resize_remap_vars(1,resize_increment)
  640. do n=1,4
  641. grid1_add_map1(num_links_old+n) = src_add(n)
  642. grid2_add_map1(num_links_old+n) = dst_add
  643. wts_map1 (:,num_links_old+n) = weights(:,n)
  644. end do
  645. case(2)
  646. num_links_old = num_links_map2
  647. num_links_map2 = num_links_old + 4
  648. if (num_links_map2 > max_links_map2)
  649. & call resize_remap_vars(2,resize_increment)
  650. do n=1,4
  651. grid1_add_map2(num_links_old+n) = dst_add
  652. grid2_add_map2(num_links_old+n) = src_add(n)
  653. wts_map2 (:,num_links_old+n) = weights(:,n)
  654. end do
  655. end select
  656. !-----------------------------------------------------------------------
  657. end subroutine store_link_bicub
  658. !***********************************************************************
  659. end module remap_bicubic
  660. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!