remap_bilinear.f 24 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781
  1. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  2. !
  3. ! this module contains necessary routines for performing an
  4. ! bilinear interpolation.
  5. !
  6. !-----------------------------------------------------------------------
  7. !
  8. ! CVS:$Id: remap_bilinear.f,v 1.6 2001/08/22 18:20:40 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_bilinear
  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_bilin
  50. !-----------------------------------------------------------------------
  51. !
  52. ! this routine computes the weights for a bilinear 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. & wgts ! bilinear weights for four corners
  70. real (kind=dbl_kind) ::
  71. & plat, plon, ! lat/lon coords of destination point
  72. & iguess, jguess, ! current guess for bilinear coordinate
  73. & thguess, phguess, ! current guess for lat/lon coordinate
  74. & deli, delj, ! corrections to i,j
  75. & dth1, dth2, dth3, ! some latitude differences
  76. & dph1, dph2, dph3, ! some longitude differences
  77. & dthp, dphp, ! difference between point and sw corner
  78. & mat1, mat2, mat3, mat4, ! matrix elements
  79. & determinant, ! matrix determinant
  80. & sum_wgts ! sum of weights for normalization
  81. !-----------------------------------------------------------------------
  82. !
  83. ! compute mappings from grid1 to grid2
  84. !
  85. !-----------------------------------------------------------------------
  86. nmap = 1
  87. if (grid1_rank /= 2) then
  88. stop 'Can not do bilinear interpolation when grid_rank /= 2'
  89. endif
  90. !***
  91. !*** loop over destination grid
  92. !***
  93. grid_loop1: do dst_add = 1, grid2_size
  94. if (.not. grid2_mask(dst_add)) cycle grid_loop1
  95. plat = grid2_center_lat(dst_add)
  96. plon = grid2_center_lon(dst_add)
  97. !***
  98. !*** find nearest square of grid points on source grid
  99. !***
  100. call grid_search_bilin(src_add, src_lats, src_lons,
  101. & plat, plon, grid1_dims,
  102. & grid1_center_lat, grid1_center_lon,
  103. & grid1_bound_box, bin_addr1, bin_addr2)
  104. !***
  105. !*** check to see if points are land points
  106. !***
  107. if (src_add(1) > 0) then
  108. do n=1,4
  109. if (.not. grid1_mask(src_add(n))) src_add(1) = 0
  110. end do
  111. endif
  112. !***
  113. !*** if point found, find local i,j coordinates for weights
  114. !***
  115. if (src_add(1) > 0) then
  116. grid2_frac(dst_add) = one
  117. !***
  118. !*** iterate to find i,j for bilinear approximation
  119. !***
  120. dth1 = src_lats(2) - src_lats(1)
  121. dth2 = src_lats(4) - src_lats(1)
  122. dth3 = src_lats(3) - src_lats(2) - dth2
  123. dph1 = src_lons(2) - src_lons(1)
  124. dph2 = src_lons(4) - src_lons(1)
  125. dph3 = src_lons(3) - src_lons(2)
  126. if (dph1 > three*pih) dph1 = dph1 - pi2
  127. if (dph2 > three*pih) dph2 = dph2 - pi2
  128. if (dph3 > three*pih) dph3 = dph3 - pi2
  129. if (dph1 < -three*pih) dph1 = dph1 + pi2
  130. if (dph2 < -three*pih) dph2 = dph2 + pi2
  131. if (dph3 < -three*pih) dph3 = dph3 + pi2
  132. dph3 = dph3 - dph2
  133. iguess = half
  134. jguess = half
  135. iter_loop1: do iter=1,max_iter
  136. dthp = plat - src_lats(1) - dth1*iguess -
  137. & dth2*jguess - dth3*iguess*jguess
  138. dphp = plon - src_lons(1)
  139. if (dphp > three*pih) dphp = dphp - pi2
  140. if (dphp < -three*pih) dphp = dphp + pi2
  141. dphp = dphp - dph1*iguess - dph2*jguess -
  142. & dph3*iguess*jguess
  143. mat1 = dth1 + dth3*jguess
  144. mat2 = dth2 + dth3*iguess
  145. mat3 = dph1 + dph3*jguess
  146. mat4 = dph2 + dph3*iguess
  147. determinant = mat1*mat4 - mat2*mat3
  148. deli = (dthp*mat4 - mat2*dphp)/determinant
  149. delj = (mat1*dphp - dthp*mat3)/determinant
  150. if (abs(deli) < converge .and.
  151. & abs(delj) < converge) exit iter_loop1
  152. iguess = iguess + deli
  153. jguess = jguess + delj
  154. end do iter_loop1
  155. if (iter <= max_iter) then
  156. !***
  157. !*** successfully found i,j - compute weights
  158. !***
  159. wgts(1) = (one-iguess)*(one-jguess)
  160. wgts(2) = iguess*(one-jguess)
  161. wgts(3) = iguess*jguess
  162. wgts(4) = (one-iguess)*jguess
  163. call store_link_bilin(dst_add, src_add, wgts, nmap)
  164. else
  165. print *,'Point coords: ',plat,plon
  166. print *,'Dest grid lats: ',src_lats
  167. print *,'Dest grid lons: ',src_lons
  168. print *,'Dest grid addresses: ',src_add
  169. print *,'Current i,j : ',iguess, jguess
  170. stop 'Iteration for i,j exceed max iteration count'
  171. endif
  172. !***
  173. !*** search for bilinear failed - use a distance-weighted
  174. !*** average instead (this is typically near the pole)
  175. !***
  176. else if (src_add(1) < 0) then
  177. src_add = abs(src_add)
  178. icount = 0
  179. do n=1,4
  180. if (grid1_mask(src_add(n))) then
  181. icount = icount + 1
  182. else
  183. src_lats(n) = zero
  184. endif
  185. end do
  186. if (icount > 0) then
  187. !*** renormalize weights
  188. sum_wgts = sum(src_lats)
  189. wgts(1) = src_lats(1)/sum_wgts
  190. wgts(2) = src_lats(2)/sum_wgts
  191. wgts(3) = src_lats(3)/sum_wgts
  192. wgts(4) = src_lats(4)/sum_wgts
  193. grid2_frac(dst_add) = one
  194. call store_link_bilin(dst_add, src_add, wgts, nmap)
  195. endif
  196. endif
  197. end do grid_loop1
  198. !-----------------------------------------------------------------------
  199. !
  200. ! compute mappings from grid2 to grid1 if necessary
  201. !
  202. !-----------------------------------------------------------------------
  203. if (num_maps > 1) then
  204. nmap = 2
  205. if (grid2_rank /= 2) then
  206. stop 'Can not do bilinear interpolation when grid_rank /= 2'
  207. endif
  208. !***
  209. !*** loop over destination grid
  210. !***
  211. grid_loop2: do dst_add = 1, grid1_size
  212. if (.not. grid1_mask(dst_add)) cycle grid_loop2
  213. plat = grid1_center_lat(dst_add)
  214. plon = grid1_center_lon(dst_add)
  215. !***
  216. !*** find nearest square of grid points on source grid
  217. !***
  218. call grid_search_bilin(src_add, src_lats, src_lons,
  219. & plat, plon, grid2_dims,
  220. & grid2_center_lat, grid2_center_lon,
  221. & grid2_bound_box, bin_addr2, bin_addr1)
  222. !***
  223. !*** check to see if points are land points
  224. !***
  225. if (src_add(1) > 0) then
  226. do n=1,4
  227. if (.not. grid2_mask(src_add(n))) src_add(1) = 0
  228. end do
  229. endif
  230. !***
  231. !*** if point found, find i,j coordinates for weights
  232. !***
  233. if (src_add(1) > 0) then
  234. grid1_frac(dst_add) = one
  235. !***
  236. !*** iterate to find i,j for bilinear approximation
  237. !***
  238. dth1 = src_lats(2) - src_lats(1)
  239. dth2 = src_lats(4) - src_lats(1)
  240. dth3 = src_lats(3) - src_lats(2) - dth2
  241. dph1 = src_lons(2) - src_lons(1)
  242. dph2 = src_lons(4) - src_lons(1)
  243. dph3 = src_lons(3) - src_lons(2)
  244. if (dph1 > pi) dph1 = dph1 - pi2
  245. if (dph2 > pi) dph2 = dph2 - pi2
  246. if (dph3 > pi) dph3 = dph3 - pi2
  247. if (dph1 < -pi) dph1 = dph1 + pi2
  248. if (dph2 < -pi) dph2 = dph2 + pi2
  249. if (dph3 < -pi) dph3 = dph3 + pi2
  250. dph3 = dph3 - dph2
  251. iguess = zero
  252. jguess = zero
  253. iter_loop2: do iter=1,max_iter
  254. dthp = plat - src_lats(1) - dth1*iguess -
  255. & dth2*jguess - dth3*iguess*jguess
  256. dphp = plon - src_lons(1)
  257. if (dphp > pi) dphp = dphp - pi2
  258. if (dphp < -pi) dphp = dphp + pi2
  259. dphp = dphp - dph1*iguess - dph2*jguess -
  260. & dph3*iguess*jguess
  261. mat1 = dth1 + dth3*jguess
  262. mat2 = dth2 + dth3*iguess
  263. mat3 = dph1 + dph3*jguess
  264. mat4 = dph2 + dph3*iguess
  265. determinant = mat1*mat4 - mat2*mat3
  266. deli = (dthp*mat4 - mat2*dphp)/determinant
  267. delj = (mat1*dphp - dthp*mat3)/determinant
  268. if (abs(deli) < converge .and.
  269. & abs(delj) < converge) exit iter_loop2
  270. iguess = iguess + deli
  271. jguess = jguess + delj
  272. end do iter_loop2
  273. if (iter <= max_iter) then
  274. !***
  275. !*** successfully found i,j - compute weights
  276. !***
  277. wgts(1) = (one-iguess)*(one-jguess)
  278. wgts(2) = iguess*(one-jguess)
  279. wgts(3) = iguess*jguess
  280. wgts(4) = (one-iguess)*jguess
  281. call store_link_bilin(dst_add, src_add, wgts, nmap)
  282. else
  283. print *,'Point coords: ',plat,plon
  284. print *,'Dest grid lats: ',src_lats
  285. print *,'Dest grid lons: ',src_lons
  286. print *,'Dest grid addresses: ',src_add
  287. print *,'Current i,j : ',iguess, jguess
  288. stop 'Iteration for i,j exceed max iteration count'
  289. endif
  290. !***
  291. !*** search for bilinear failed - us a distance-weighted
  292. !*** average instead
  293. !***
  294. else if (src_add(1) < 0) then
  295. src_add = abs(src_add)
  296. icount = 0
  297. do n=1,4
  298. if (grid2_mask(src_add(n))) then
  299. icount = icount + 1
  300. else
  301. src_lats(n) = zero
  302. endif
  303. end do
  304. if (icount > 0) then
  305. !*** renormalize weights
  306. sum_wgts = sum(src_lats)
  307. wgts(1) = src_lats(1)/sum_wgts
  308. wgts(2) = src_lats(2)/sum_wgts
  309. wgts(3) = src_lats(3)/sum_wgts
  310. wgts(4) = src_lats(4)/sum_wgts
  311. grid1_frac(dst_add) = one
  312. call store_link_bilin(dst_add, src_add, wgts, nmap)
  313. endif
  314. endif
  315. end do grid_loop2
  316. endif ! nmap=2
  317. !-----------------------------------------------------------------------
  318. end subroutine remap_bilin
  319. !***********************************************************************
  320. subroutine grid_search_bilin(src_add, src_lats, src_lons,
  321. & plat, plon, src_grid_dims,
  322. & src_center_lat, src_center_lon,
  323. & src_grid_bound_box,
  324. & src_bin_add, dst_bin_add)
  325. !-----------------------------------------------------------------------
  326. !
  327. ! this routine finds the location of the search point plat, plon
  328. ! in the source grid and returns the corners needed for a bilinear
  329. ! interpolation.
  330. !
  331. !-----------------------------------------------------------------------
  332. !-----------------------------------------------------------------------
  333. !
  334. ! output variables
  335. !
  336. !-----------------------------------------------------------------------
  337. integer (kind=int_kind), dimension(4), intent(out) ::
  338. & src_add ! address of each corner point enclosing P
  339. real (kind=dbl_kind), dimension(4), intent(out) ::
  340. & src_lats, ! latitudes of the four corner points
  341. & src_lons ! longitudes of the four corner points
  342. !-----------------------------------------------------------------------
  343. !
  344. ! input variables
  345. !
  346. !-----------------------------------------------------------------------
  347. real (kind=dbl_kind), intent(in) ::
  348. & plat, ! latitude of the search point
  349. & plon ! longitude of the search point
  350. integer (kind=int_kind), dimension(2), intent(in) ::
  351. & src_grid_dims ! size of each src grid dimension
  352. real (kind=dbl_kind), dimension(:), intent(in) ::
  353. & src_center_lat, ! latitude of each src grid center
  354. & src_center_lon ! longitude of each src grid center
  355. real (kind=dbl_kind), dimension(:,:), intent(in) ::
  356. & src_grid_bound_box ! bound box for source grid
  357. integer (kind=int_kind), dimension(:,:), intent(in) ::
  358. & src_bin_add, ! latitude bins for restricting
  359. & dst_bin_add ! searches
  360. !-----------------------------------------------------------------------
  361. !
  362. ! local variables
  363. !
  364. !-----------------------------------------------------------------------
  365. integer (kind=int_kind) :: n, next_n, srch_add, ! dummy indices
  366. & nx, ny, ! dimensions of src grid
  367. & min_add, max_add, ! addresses for restricting search
  368. & i, j, jp1, ip1, n_add, e_add, ne_add ! addresses
  369. real (kind=dbl_kind) :: ! vectors for cross-product check
  370. & vec1_lat, vec1_lon,
  371. & vec2_lat, vec2_lon, cross_product, cross_product_last,
  372. & coslat_dst, sinlat_dst, coslon_dst, sinlon_dst,
  373. & dist_min, distance ! for computing dist-weighted avg
  374. !-----------------------------------------------------------------------
  375. !
  376. ! restrict search first using bins
  377. !
  378. !-----------------------------------------------------------------------
  379. src_add = 0
  380. min_add = size(src_center_lat)
  381. max_add = 1
  382. do n=1,num_srch_bins
  383. if (plat >= bin_lats(1,n) .and. plat <= bin_lats(2,n) .and.
  384. & plon >= bin_lons(1,n) .and. plon <= bin_lons(2,n)) then
  385. min_add = min(min_add, src_bin_add(1,n))
  386. max_add = max(max_add, src_bin_add(2,n))
  387. endif
  388. end do
  389. !-----------------------------------------------------------------------
  390. !
  391. ! now perform a more detailed search
  392. !
  393. !-----------------------------------------------------------------------
  394. nx = src_grid_dims(1)
  395. ny = src_grid_dims(2)
  396. srch_loop: do srch_add = min_add,max_add
  397. !*** first check bounding box
  398. if (plat <= src_grid_bound_box(2,srch_add) .and.
  399. & plat >= src_grid_bound_box(1,srch_add) .and.
  400. & plon <= src_grid_bound_box(4,srch_add) .and.
  401. & plon >= src_grid_bound_box(3,srch_add)) then
  402. !***
  403. !*** we are within bounding box so get really serious
  404. !***
  405. !*** determine neighbor addresses
  406. j = (srch_add - 1)/nx +1
  407. i = srch_add - (j-1)*nx
  408. if (i < nx) then
  409. ip1 = i + 1
  410. else
  411. ip1 = 1
  412. endif
  413. if (j < ny) then
  414. jp1 = j+1
  415. else
  416. jp1 = 1
  417. endif
  418. n_add = (jp1 - 1)*nx + i
  419. e_add = (j - 1)*nx + ip1
  420. ne_add = (jp1 - 1)*nx + ip1
  421. src_lats(1) = src_center_lat(srch_add)
  422. src_lats(2) = src_center_lat(e_add)
  423. src_lats(3) = src_center_lat(ne_add)
  424. src_lats(4) = src_center_lat(n_add)
  425. src_lons(1) = src_center_lon(srch_add)
  426. src_lons(2) = src_center_lon(e_add)
  427. src_lons(3) = src_center_lon(ne_add)
  428. src_lons(4) = src_center_lon(n_add)
  429. !***
  430. !*** for consistency, we must make sure all lons are in
  431. !*** same 2pi interval
  432. !***
  433. vec1_lon = src_lons(1) - plon
  434. if (vec1_lon > pi) then
  435. src_lons(1) = src_lons(1) - pi2
  436. else if (vec1_lon < -pi) then
  437. src_lons(1) = src_lons(1) + pi2
  438. endif
  439. do n=2,4
  440. vec1_lon = src_lons(n) - src_lons(1)
  441. if (vec1_lon > pi) then
  442. src_lons(n) = src_lons(n) - pi2
  443. else if (vec1_lon < -pi) then
  444. src_lons(n) = src_lons(n) + pi2
  445. endif
  446. end do
  447. corner_loop: do n=1,4
  448. next_n = MOD(n,4) + 1
  449. !***
  450. !*** here we take the cross product of the vector making
  451. !*** up each box side with the vector formed by the vertex
  452. !*** and search point. if all the cross products are
  453. !*** positive, the point is contained in the box.
  454. !***
  455. vec1_lat = src_lats(next_n) - src_lats(n)
  456. vec1_lon = src_lons(next_n) - src_lons(n)
  457. vec2_lat = plat - src_lats(n)
  458. vec2_lon = plon - src_lons(n)
  459. !***
  460. !*** check for 0,2pi crossings
  461. !***
  462. if (vec1_lon > three*pih) then
  463. vec1_lon = vec1_lon - pi2
  464. else if (vec1_lon < -three*pih) then
  465. vec1_lon = vec1_lon + pi2
  466. endif
  467. if (vec2_lon > three*pih) then
  468. vec2_lon = vec2_lon - pi2
  469. else if (vec2_lon < -three*pih) then
  470. vec2_lon = vec2_lon + pi2
  471. endif
  472. cross_product = vec1_lon*vec2_lat - vec2_lon*vec1_lat
  473. !***
  474. !*** if cross product is less than zero, this cell
  475. !*** doesn't work
  476. !***
  477. if (n == 1) cross_product_last = cross_product
  478. if (cross_product*cross_product_last < zero)
  479. & exit corner_loop
  480. cross_product_last = cross_product
  481. end do corner_loop
  482. !***
  483. !*** if cross products all same sign, we found the location
  484. !***
  485. if (n > 4) then
  486. src_add(1) = srch_add
  487. src_add(2) = e_add
  488. src_add(3) = ne_add
  489. src_add(4) = n_add
  490. return
  491. endif
  492. !***
  493. !*** otherwise move on to next cell
  494. !***
  495. endif !bounding box check
  496. end do srch_loop
  497. !***
  498. !*** if no cell found, point is likely either in a box that
  499. !*** straddles either pole or is outside the grid. fall back
  500. !*** to a distance-weighted average of the four closest
  501. !*** points. go ahead and compute weights here, but store
  502. !*** in src_lats and return -add to prevent the parent
  503. !*** routine from computing bilinear weights
  504. !***
  505. !print *,'Could not find location for ',plat,plon
  506. !print *,'Using nearest-neighbor average for this point'
  507. coslat_dst = cos(plat)
  508. sinlat_dst = sin(plat)
  509. coslon_dst = cos(plon)
  510. sinlon_dst = sin(plon)
  511. dist_min = bignum
  512. src_lats = bignum
  513. do srch_add = min_add,max_add
  514. distance = acos(coslat_dst*cos(src_center_lat(srch_add))*
  515. & (coslon_dst*cos(src_center_lon(srch_add)) +
  516. & sinlon_dst*sin(src_center_lon(srch_add)))+
  517. & sinlat_dst*sin(src_center_lat(srch_add)))
  518. if (distance < dist_min) then
  519. sort_loop: do n=1,4
  520. if (distance < src_lats(n)) then
  521. do i=4,n+1,-1
  522. src_add (i) = src_add (i-1)
  523. src_lats(i) = src_lats(i-1)
  524. end do
  525. src_add (n) = -srch_add
  526. src_lats(n) = distance
  527. dist_min = src_lats(4)
  528. exit sort_loop
  529. endif
  530. end do sort_loop
  531. endif
  532. end do
  533. src_lons = one/(src_lats + tiny)
  534. distance = sum(src_lons)
  535. src_lats = src_lons/distance
  536. !-----------------------------------------------------------------------
  537. end subroutine grid_search_bilin
  538. !***********************************************************************
  539. subroutine store_link_bilin(dst_add, src_add, weights, nmap)
  540. !-----------------------------------------------------------------------
  541. !
  542. ! this routine stores the address and weight for four links
  543. ! associated with one destination point in the appropriate address
  544. ! and weight arrays and resizes those arrays if necessary.
  545. !
  546. !-----------------------------------------------------------------------
  547. !-----------------------------------------------------------------------
  548. !
  549. ! input variables
  550. !
  551. !-----------------------------------------------------------------------
  552. integer (kind=int_kind), intent(in) ::
  553. & dst_add, ! address on destination grid
  554. & nmap ! identifies which direction for mapping
  555. integer (kind=int_kind), dimension(4), intent(in) ::
  556. & src_add ! addresses on source grid
  557. real (kind=dbl_kind), dimension(4), intent(in) ::
  558. & weights ! array of remapping weights for these links
  559. !-----------------------------------------------------------------------
  560. !
  561. ! local variables
  562. !
  563. !-----------------------------------------------------------------------
  564. integer (kind=int_kind) :: n, ! dummy index
  565. & num_links_old ! placeholder for old link number
  566. !-----------------------------------------------------------------------
  567. !
  568. ! increment number of links and check to see if remap arrays need
  569. ! to be increased to accomodate the new link. then store the
  570. ! link.
  571. !
  572. !-----------------------------------------------------------------------
  573. select case (nmap)
  574. case(1)
  575. num_links_old = num_links_map1
  576. num_links_map1 = num_links_old + 4
  577. if (num_links_map1 > max_links_map1)
  578. & call resize_remap_vars(1,resize_increment)
  579. do n=1,4
  580. grid1_add_map1(num_links_old+n) = src_add(n)
  581. grid2_add_map1(num_links_old+n) = dst_add
  582. wts_map1 (1,num_links_old+n) = weights(n)
  583. end do
  584. case(2)
  585. num_links_old = num_links_map2
  586. num_links_map2 = num_links_old + 4
  587. if (num_links_map2 > max_links_map2)
  588. & call resize_remap_vars(2,resize_increment)
  589. do n=1,4
  590. grid1_add_map2(num_links_old+n) = dst_add
  591. grid2_add_map2(num_links_old+n) = src_add(n)
  592. wts_map2 (1,num_links_old+n) = weights(n)
  593. end do
  594. end select
  595. !-----------------------------------------------------------------------
  596. end subroutine store_link_bilin
  597. !***********************************************************************
  598. end module remap_bilinear
  599. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!