remap_conserv.f 73 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396139713981399140014011402140314041405140614071408140914101411141214131414141514161417141814191420142114221423142414251426142714281429143014311432143314341435143614371438143914401441144214431444144514461447144814491450145114521453145414551456145714581459146014611462146314641465146614671468146914701471147214731474147514761477147814791480148114821483148414851486148714881489149014911492149314941495149614971498149915001501150215031504150515061507150815091510151115121513151415151516151715181519152015211522152315241525152615271528152915301531153215331534153515361537153815391540154115421543154415451546154715481549155015511552155315541555155615571558155915601561156215631564156515661567156815691570157115721573157415751576157715781579158015811582158315841585158615871588158915901591159215931594159515961597159815991600160116021603160416051606160716081609161016111612161316141615161616171618161916201621162216231624162516261627162816291630163116321633163416351636163716381639164016411642164316441645164616471648164916501651165216531654165516561657165816591660166116621663166416651666166716681669167016711672167316741675167616771678167916801681168216831684168516861687168816891690169116921693169416951696169716981699170017011702170317041705170617071708170917101711171217131714171517161717171817191720172117221723172417251726172717281729173017311732173317341735173617371738173917401741174217431744174517461747174817491750175117521753175417551756175717581759176017611762176317641765176617671768176917701771177217731774177517761777177817791780178117821783178417851786178717881789179017911792179317941795179617971798179918001801180218031804180518061807180818091810181118121813181418151816181718181819182018211822182318241825182618271828182918301831183218331834183518361837183818391840184118421843184418451846184718481849185018511852185318541855185618571858185918601861186218631864186518661867186818691870187118721873187418751876187718781879188018811882188318841885188618871888188918901891189218931894189518961897189818991900190119021903190419051906190719081909191019111912191319141915191619171918191919201921192219231924192519261927192819291930193119321933193419351936193719381939194019411942194319441945194619471948194919501951195219531954195519561957195819591960196119621963196419651966196719681969197019711972197319741975197619771978197919801981198219831984198519861987198819891990199119921993199419951996199719981999200020012002200320042005200620072008200920102011201220132014201520162017201820192020202120222023202420252026202720282029203020312032203320342035203620372038203920402041204220432044204520462047204820492050205120522053205420552056205720582059206020612062206320642065206620672068206920702071207220732074207520762077207820792080208120822083208420852086208720882089209020912092209320942095209620972098209921002101210221032104210521062107210821092110211121122113211421152116211721182119212021212122212321242125212621272128212921302131213221332134213521362137213821392140214121422143214421452146214721482149215021512152215321542155215621572158215921602161216221632164216521662167216821692170217121722173217421752176217721782179218021812182218321842185218621872188218921902191219221932194219521962197
  1. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  2. !
  3. ! this module contains necessary routines for computing addresses
  4. ! and weights for a conservative interpolation between any two
  5. ! grids on a sphere. the weights are computed by performing line
  6. ! integrals around all overlap regions of the two grids. see
  7. ! Dukowicz and Kodis, SIAM J. Sci. Stat. Comput. 8, 305 (1987) and
  8. ! Jones, P.W. Monthly Weather Review (submitted).
  9. !
  10. !-----------------------------------------------------------------------
  11. !
  12. ! CVS:$Id: remap_conserv.f,v 1.10 2001/08/21 21:05:13 pwjones Exp $
  13. !
  14. ! Copyright (c) 1997, 1998 the Regents of the University of
  15. ! California.
  16. !
  17. ! This software and ancillary information (herein called software)
  18. ! called SCRIP is made available under the terms described here.
  19. ! The software has been approved for release with associated
  20. ! LA-CC Number 98-45.
  21. !
  22. ! Unless otherwise indicated, this software has been authored
  23. ! by an employee or employees of the University of California,
  24. ! operator of the Los Alamos National Laboratory under Contract
  25. ! No. W-7405-ENG-36 with the U.S. Department of Energy. The U.S.
  26. ! Government has rights to use, reproduce, and distribute this
  27. ! software. The public may copy and use this software without
  28. ! charge, provided that this Notice and any statement of authorship
  29. ! are reproduced on all copies. Neither the Government nor the
  30. ! University makes any warranty, express or implied, or assumes
  31. ! any liability or responsibility for the use of this software.
  32. !
  33. ! If software is modified to produce derivative works, such modified
  34. ! software should be clearly marked, so as not to confuse it with
  35. ! the version available from Los Alamos National Laboratory.
  36. !
  37. !***********************************************************************
  38. module remap_conservative
  39. !-----------------------------------------------------------------------
  40. use kinds_mod ! defines common data types
  41. use constants ! defines common constants
  42. use timers ! module for timing
  43. use grids ! module containing grid information
  44. use remap_vars ! module containing remap information
  45. implicit none
  46. !-----------------------------------------------------------------------
  47. !
  48. ! module variables
  49. !
  50. !-----------------------------------------------------------------------
  51. integer (kind=int_kind), save ::
  52. & num_srch_cells ! num cells in restricted search arrays
  53. integer (kind=int_kind), dimension(:), allocatable, save ::
  54. & srch_add ! global address of cells in srch arrays
  55. real (kind=dbl_kind), parameter ::
  56. & north_thresh = 1.45_dbl_kind, ! threshold for coord transf.
  57. & south_thresh =-2.00_dbl_kind ! threshold for coord transf.
  58. real (kind=dbl_kind), dimension(:,:), allocatable, save ::
  59. & srch_corner_lat, ! lat of each corner of srch cells
  60. & srch_corner_lon ! lon of each corner of srch cells
  61. !***********************************************************************
  62. contains
  63. !***********************************************************************
  64. subroutine remap_conserv
  65. !-----------------------------------------------------------------------
  66. !
  67. ! this routine traces the perimeters of every grid cell on each
  68. ! grid checking for intersections with the other grid and computing
  69. ! line integrals for each subsegment.
  70. !
  71. !-----------------------------------------------------------------------
  72. !-----------------------------------------------------------------------
  73. !
  74. ! local variables
  75. !
  76. !-----------------------------------------------------------------------
  77. integer (kind=int_kind), parameter ::
  78. & max_subseg = 10000 ! max number of subsegments per segment
  79. ! to prevent infinite loop
  80. integer (kind=int_kind) ::
  81. & grid1_add, ! current linear address for grid1 cell
  82. & grid2_add, ! current linear address for grid2 cell
  83. & min_add, ! addresses for restricting search of
  84. & max_add, ! destination grid
  85. & n, nwgt, ! generic counters
  86. & corner, ! corner of cell that segment starts from
  87. & next_corn, ! corner of cell that segment ends on
  88. & num_subseg ! number of subsegments
  89. logical (kind=log_kind) ::
  90. & lcoinc, ! flag for coincident segments
  91. & lrevers, ! flag for reversing direction of segment
  92. & lbegin ! flag for first integration of a segment
  93. logical (kind=log_kind), dimension(:), allocatable ::
  94. & srch_mask ! mask for restricting searches
  95. real (kind=dbl_kind) ::
  96. & intrsct_lat, intrsct_lon, ! lat/lon of next intersect
  97. & beglat, endlat, beglon, endlon, ! endpoints of current seg.
  98. & norm_factor ! factor for normalizing wts
  99. real (kind=dbl_kind), dimension(:), allocatable ::
  100. & grid2_centroid_lat, grid2_centroid_lon, ! centroid coords
  101. & grid1_centroid_lat, grid1_centroid_lon ! on each grid
  102. real (kind=dbl_kind), dimension(2) :: begseg ! begin lat/lon for
  103. ! full segment
  104. real (kind=dbl_kind), dimension(6) :: weights ! local wgt array
  105. !-----------------------------------------------------------------------
  106. !
  107. ! initialize centroid arrays
  108. !
  109. !-----------------------------------------------------------------------
  110. allocate( grid1_centroid_lat(grid1_size),
  111. & grid1_centroid_lon(grid1_size),
  112. & grid2_centroid_lat(grid2_size),
  113. & grid2_centroid_lon(grid2_size))
  114. grid1_centroid_lat = zero
  115. grid1_centroid_lon = zero
  116. grid2_centroid_lat = zero
  117. grid2_centroid_lon = zero
  118. !-----------------------------------------------------------------------
  119. !
  120. ! integrate around each cell on grid1
  121. !
  122. !-----------------------------------------------------------------------
  123. allocate(srch_mask(grid2_size))
  124. print *,'grid1 sweep '
  125. do grid1_add = 1,grid1_size
  126. !***
  127. !*** restrict searches first using search bins
  128. !***
  129. call timer_start(1)
  130. min_add = grid2_size
  131. max_add = 1
  132. do n=1,num_srch_bins
  133. if (grid1_add >= bin_addr1(1,n) .and.
  134. & grid1_add <= bin_addr1(2,n)) then
  135. min_add = min(min_add, bin_addr2(1,n))
  136. max_add = max(max_add, bin_addr2(2,n))
  137. endif
  138. end do
  139. !***
  140. !*** further restrict searches using bounding boxes
  141. !***
  142. num_srch_cells = 0
  143. do grid2_add = min_add,max_add
  144. srch_mask(grid2_add) = (grid2_bound_box(1,grid2_add) <=
  145. & grid1_bound_box(2,grid1_add)) .and.
  146. & (grid2_bound_box(2,grid2_add) >=
  147. & grid1_bound_box(1,grid1_add)) .and.
  148. & (grid2_bound_box(3,grid2_add) <=
  149. & grid1_bound_box(4,grid1_add)) .and.
  150. & (grid2_bound_box(4,grid2_add) >=
  151. & grid1_bound_box(3,grid1_add))
  152. if (srch_mask(grid2_add)) num_srch_cells = num_srch_cells+1
  153. end do
  154. !***
  155. !*** create search arrays
  156. !***
  157. allocate(srch_add(num_srch_cells),
  158. & srch_corner_lat(grid2_corners,num_srch_cells),
  159. & srch_corner_lon(grid2_corners,num_srch_cells))
  160. n = 0
  161. gather1: do grid2_add = min_add,max_add
  162. if (srch_mask(grid2_add)) then
  163. n = n+1
  164. srch_add(n) = grid2_add
  165. srch_corner_lat(:,n) = grid2_corner_lat(:,grid2_add)
  166. srch_corner_lon(:,n) = grid2_corner_lon(:,grid2_add)
  167. endif
  168. end do gather1
  169. call timer_stop(1)
  170. !***
  171. !*** integrate around this cell
  172. !***
  173. do corner = 1,grid1_corners
  174. next_corn = mod(corner,grid1_corners) + 1
  175. !***
  176. !*** define endpoints of the current segment
  177. !***
  178. beglat = grid1_corner_lat(corner,grid1_add)
  179. beglon = grid1_corner_lon(corner,grid1_add)
  180. endlat = grid1_corner_lat(next_corn,grid1_add)
  181. endlon = grid1_corner_lon(next_corn,grid1_add)
  182. lrevers = .false.
  183. !***
  184. !*** to ensure exact path taken during both
  185. !*** sweeps, always integrate segments in the same
  186. !*** direction (SW to NE).
  187. !***
  188. if ((endlat < beglat) .or.
  189. & (endlat == beglat .and. endlon < beglon)) then
  190. beglat = grid1_corner_lat(next_corn,grid1_add)
  191. beglon = grid1_corner_lon(next_corn,grid1_add)
  192. endlat = grid1_corner_lat(corner,grid1_add)
  193. endlon = grid1_corner_lon(corner,grid1_add)
  194. lrevers = .true.
  195. endif
  196. begseg(1) = beglat
  197. begseg(2) = beglon
  198. lbegin = .true.
  199. num_subseg = 0
  200. !***
  201. !*** if this is a constant-longitude segment, skip the rest
  202. !*** since the line integral contribution will be zero.
  203. !***
  204. if (endlon /= beglon) then
  205. !***
  206. !*** integrate along this segment, detecting intersections
  207. !*** and computing the line integral for each sub-segment
  208. !***
  209. do while (beglat /= endlat .or. beglon /= endlon)
  210. !***
  211. !*** prevent infinite loops if integration gets stuck
  212. !*** near cell or threshold boundary
  213. !***
  214. num_subseg = num_subseg + 1
  215. if (num_subseg > max_subseg) then
  216. stop 'integration stalled: num_subseg exceeded limit'
  217. endif
  218. !***
  219. !*** find next intersection of this segment with a grid
  220. !*** line on grid 2.
  221. !***
  222. call timer_start(2)
  223. call intersection(grid2_add,intrsct_lat,intrsct_lon,lcoinc,
  224. & beglat, beglon, endlat, endlon, begseg,
  225. & lbegin, lrevers)
  226. call timer_stop(2)
  227. lbegin = .false.
  228. !***
  229. !*** compute line integral for this subsegment.
  230. !***
  231. call timer_start(3)
  232. if (grid2_add /= 0) then
  233. call line_integral(weights, num_wts,
  234. & beglon, intrsct_lon, beglat, intrsct_lat,
  235. & grid1_center_lat(grid1_add),
  236. & grid1_center_lon(grid1_add),
  237. & grid2_center_lat(grid2_add),
  238. & grid2_center_lon(grid2_add))
  239. else
  240. call line_integral(weights, num_wts,
  241. & beglon, intrsct_lon, beglat, intrsct_lat,
  242. & grid1_center_lat(grid1_add),
  243. & grid1_center_lon(grid1_add),
  244. & grid1_center_lat(grid1_add),
  245. & grid1_center_lon(grid1_add))
  246. endif
  247. call timer_stop(3)
  248. !***
  249. !*** if integrating in reverse order, change
  250. !*** sign of weights
  251. !***
  252. if (lrevers) then
  253. weights = -weights
  254. endif
  255. !***
  256. !*** store the appropriate addresses and weights.
  257. !*** also add contributions to cell areas and centroids.
  258. !***
  259. !if (grid1_add == 119247) then
  260. ! print *,grid1_add,grid2_add,corner,weights(1)
  261. ! print *,grid1_corner_lat(:,grid1_add)
  262. ! print *,grid1_corner_lon(:,grid1_add)
  263. ! print *,grid2_corner_lat(:,grid2_add)
  264. ! print *,grid2_corner_lon(:,grid2_add)
  265. ! print *,beglat,beglon,intrsct_lat,intrsct_lon
  266. !endif
  267. if (grid2_add /= 0) then
  268. if (grid1_mask(grid1_add)) then
  269. call timer_start(4)
  270. call store_link_cnsrv(grid1_add, grid2_add, weights)
  271. call timer_stop(4)
  272. grid1_frac(grid1_add) = grid1_frac(grid1_add) +
  273. & weights(1)
  274. grid2_frac(grid2_add) = grid2_frac(grid2_add) +
  275. & weights(num_wts+1)
  276. endif
  277. endif
  278. grid1_area(grid1_add) = grid1_area(grid1_add) + weights(1)
  279. grid1_centroid_lat(grid1_add) =
  280. & grid1_centroid_lat(grid1_add) + weights(2)
  281. grid1_centroid_lon(grid1_add) =
  282. & grid1_centroid_lon(grid1_add) + weights(3)
  283. !***
  284. !*** reset beglat and beglon for next subsegment.
  285. !***
  286. beglat = intrsct_lat
  287. beglon = intrsct_lon
  288. end do
  289. endif
  290. !***
  291. !*** end of segment
  292. !***
  293. end do
  294. !***
  295. !*** finished with this cell: deallocate search array and
  296. !*** start on next cell
  297. deallocate(srch_add, srch_corner_lat, srch_corner_lon)
  298. end do
  299. deallocate(srch_mask)
  300. !-----------------------------------------------------------------------
  301. !
  302. ! integrate around each cell on grid2
  303. !
  304. !-----------------------------------------------------------------------
  305. allocate(srch_mask(grid1_size))
  306. print *,'grid2 sweep '
  307. do grid2_add = 1,grid2_size
  308. !***
  309. !*** restrict searches first using search bins
  310. !***
  311. call timer_start(5)
  312. min_add = grid1_size
  313. max_add = 1
  314. do n=1,num_srch_bins
  315. if (grid2_add >= bin_addr2(1,n) .and.
  316. & grid2_add <= bin_addr2(2,n)) then
  317. min_add = min(min_add, bin_addr1(1,n))
  318. max_add = max(max_add, bin_addr1(2,n))
  319. endif
  320. end do
  321. !***
  322. !*** further restrict searches using bounding boxes
  323. !***
  324. num_srch_cells = 0
  325. do grid1_add = min_add, max_add
  326. srch_mask(grid1_add) = (grid1_bound_box(1,grid1_add) <=
  327. & grid2_bound_box(2,grid2_add)) .and.
  328. & (grid1_bound_box(2,grid1_add) >=
  329. & grid2_bound_box(1,grid2_add)) .and.
  330. & (grid1_bound_box(3,grid1_add) <=
  331. & grid2_bound_box(4,grid2_add)) .and.
  332. & (grid1_bound_box(4,grid1_add) >=
  333. & grid2_bound_box(3,grid2_add))
  334. if (srch_mask(grid1_add)) num_srch_cells = num_srch_cells+1
  335. end do
  336. allocate(srch_add(num_srch_cells),
  337. & srch_corner_lat(grid1_corners,num_srch_cells),
  338. & srch_corner_lon(grid1_corners,num_srch_cells))
  339. n = 0
  340. gather2: do grid1_add = min_add,max_add
  341. if (srch_mask(grid1_add)) then
  342. n = n+1
  343. srch_add(n) = grid1_add
  344. srch_corner_lat(:,n) = grid1_corner_lat(:,grid1_add)
  345. srch_corner_lon(:,n) = grid1_corner_lon(:,grid1_add)
  346. endif
  347. end do gather2
  348. call timer_stop(5)
  349. !***
  350. !*** integrate around this cell
  351. !***
  352. do corner = 1,grid2_corners
  353. next_corn = mod(corner,grid2_corners) + 1
  354. beglat = grid2_corner_lat(corner,grid2_add)
  355. beglon = grid2_corner_lon(corner,grid2_add)
  356. endlat = grid2_corner_lat(next_corn,grid2_add)
  357. endlon = grid2_corner_lon(next_corn,grid2_add)
  358. lrevers = .false.
  359. !***
  360. !*** to ensure exact path taken during both
  361. !*** sweeps, always integrate in the same direction
  362. !***
  363. if ((endlat < beglat) .or.
  364. & (endlat == beglat .and. endlon < beglon)) then
  365. beglat = grid2_corner_lat(next_corn,grid2_add)
  366. beglon = grid2_corner_lon(next_corn,grid2_add)
  367. endlat = grid2_corner_lat(corner,grid2_add)
  368. endlon = grid2_corner_lon(corner,grid2_add)
  369. lrevers = .true.
  370. endif
  371. begseg(1) = beglat
  372. begseg(2) = beglon
  373. lbegin = .true.
  374. !***
  375. !*** if this is a constant-longitude segment, skip the rest
  376. !*** since the line integral contribution will be zero.
  377. !***
  378. if (endlon /= beglon) then
  379. num_subseg = 0
  380. !***
  381. !*** integrate along this segment, detecting intersections
  382. !*** and computing the line integral for each sub-segment
  383. !***
  384. do while (beglat /= endlat .or. beglon /= endlon)
  385. !***
  386. !*** prevent infinite loops if integration gets stuck
  387. !*** near cell or threshold boundary
  388. !***
  389. num_subseg = num_subseg + 1
  390. if (num_subseg > max_subseg) then
  391. stop 'integration stalled: num_subseg exceeded limit'
  392. endif
  393. !***
  394. !*** find next intersection of this segment with a line
  395. !*** on grid 2.
  396. !***
  397. call timer_start(6)
  398. call intersection(grid1_add,intrsct_lat,intrsct_lon,lcoinc,
  399. & beglat, beglon, endlat, endlon, begseg,
  400. & lbegin, lrevers)
  401. call timer_stop(6)
  402. lbegin = .false.
  403. !***
  404. !*** compute line integral for this subsegment.
  405. !***
  406. call timer_start(7)
  407. if (grid1_add /= 0) then
  408. call line_integral(weights, num_wts,
  409. & beglon, intrsct_lon, beglat, intrsct_lat,
  410. & grid1_center_lat(grid1_add),
  411. & grid1_center_lon(grid1_add),
  412. & grid2_center_lat(grid2_add),
  413. & grid2_center_lon(grid2_add))
  414. else
  415. call line_integral(weights, num_wts,
  416. & beglon, intrsct_lon, beglat, intrsct_lat,
  417. & grid2_center_lat(grid2_add),
  418. & grid2_center_lon(grid2_add),
  419. & grid2_center_lat(grid2_add),
  420. & grid2_center_lon(grid2_add))
  421. endif
  422. call timer_stop(7)
  423. if (lrevers) then
  424. weights = -weights
  425. endif
  426. !***
  427. !*** store the appropriate addresses and weights.
  428. !*** also add contributions to cell areas and centroids.
  429. !*** if there is a coincidence, do not store weights
  430. !*** because they have been captured in the previous loop.
  431. !*** the grid1 mask is the master mask
  432. !***
  433. !if (grid1_add == 119247) then
  434. ! print *,grid1_add,grid2_add,corner,weights(1)
  435. ! print *,grid1_corner_lat(:,grid1_add)
  436. ! print *,grid1_corner_lon(:,grid1_add)
  437. ! print *,grid2_corner_lat(:,grid2_add)
  438. ! print *,grid2_corner_lon(:,grid2_add)
  439. ! print *,beglat,beglon,intrsct_lat,intrsct_lon
  440. !endif
  441. if (.not. lcoinc .and. grid1_add /= 0) then
  442. if (grid1_mask(grid1_add)) then
  443. call timer_start(8)
  444. call store_link_cnsrv(grid1_add, grid2_add, weights)
  445. call timer_stop(8)
  446. grid1_frac(grid1_add) = grid1_frac(grid1_add) +
  447. & weights(1)
  448. grid2_frac(grid2_add) = grid2_frac(grid2_add) +
  449. & weights(num_wts+1)
  450. endif
  451. endif
  452. grid2_area(grid2_add) = grid2_area(grid2_add) +
  453. & weights(num_wts+1)
  454. grid2_centroid_lat(grid2_add) =
  455. & grid2_centroid_lat(grid2_add) + weights(num_wts+2)
  456. grid2_centroid_lon(grid2_add) =
  457. & grid2_centroid_lon(grid2_add) + weights(num_wts+3)
  458. !***
  459. !*** reset beglat and beglon for next subsegment.
  460. !***
  461. beglat = intrsct_lat
  462. beglon = intrsct_lon
  463. end do
  464. endif
  465. !***
  466. !*** end of segment
  467. !***
  468. end do
  469. !***
  470. !*** finished with this cell: deallocate search array and
  471. !*** start on next cell
  472. deallocate(srch_add, srch_corner_lat, srch_corner_lon)
  473. end do
  474. deallocate(srch_mask)
  475. !-----------------------------------------------------------------------
  476. !
  477. ! correct for situations where N/S pole not explicitly included in
  478. ! grid (i.e. as a grid corner point). if pole is missing from only
  479. ! one grid, need to correct only the area and centroid of that
  480. ! grid. if missing from both, do complete weight calculation.
  481. !
  482. !-----------------------------------------------------------------------
  483. !*** North Pole
  484. weights(1) = pi2
  485. weights(2) = pi*pi
  486. weights(3) = zero
  487. weights(4) = pi2
  488. weights(5) = pi*pi
  489. weights(6) = zero
  490. grid1_add = 0
  491. pole_loop1: do n=1,grid1_size
  492. if (grid1_area(n) < -three*pih .and.
  493. & grid1_center_lat(n) > zero) then
  494. grid1_add = n
  495. exit pole_loop1
  496. endif
  497. end do pole_loop1
  498. grid2_add = 0
  499. pole_loop2: do n=1,grid2_size
  500. if (grid2_area(n) < -three*pih .and.
  501. & grid2_center_lat(n) > zero) then
  502. grid2_add = n
  503. exit pole_loop2
  504. endif
  505. end do pole_loop2
  506. if (grid1_add /=0) then
  507. grid1_area(grid1_add) = grid1_area(grid1_add) + weights(1)
  508. grid1_centroid_lat(grid1_add) =
  509. & grid1_centroid_lat(grid1_add) + weights(2)
  510. grid1_centroid_lon(grid1_add) =
  511. & grid1_centroid_lon(grid1_add) + weights(3)
  512. endif
  513. if (grid2_add /=0) then
  514. grid2_area(grid2_add) = grid2_area(grid2_add) +
  515. & weights(num_wts+1)
  516. grid2_centroid_lat(grid2_add) =
  517. & grid2_centroid_lat(grid2_add) + weights(num_wts+2)
  518. grid2_centroid_lon(grid2_add) =
  519. & grid2_centroid_lon(grid2_add) + weights(num_wts+3)
  520. endif
  521. if (grid1_add /= 0 .and. grid2_add /=0) then
  522. call store_link_cnsrv(grid1_add, grid2_add, weights)
  523. grid1_frac(grid1_add) = grid1_frac(grid1_add) +
  524. & weights(1)
  525. grid2_frac(grid2_add) = grid2_frac(grid2_add) +
  526. & weights(num_wts+1)
  527. endif
  528. !*** South Pole
  529. weights(1) = pi2
  530. weights(2) = -pi*pi
  531. weights(3) = zero
  532. weights(4) = pi2
  533. weights(5) = -pi*pi
  534. weights(6) = zero
  535. grid1_add = 0
  536. pole_loop3: do n=1,grid1_size
  537. if (grid1_area(n) < -three*pih .and.
  538. & grid1_center_lat(n) < zero) then
  539. grid1_add = n
  540. exit pole_loop3
  541. endif
  542. end do pole_loop3
  543. grid2_add = 0
  544. pole_loop4: do n=1,grid2_size
  545. if (grid2_area(n) < -three*pih .and.
  546. & grid2_center_lat(n) < zero) then
  547. grid2_add = n
  548. exit pole_loop4
  549. endif
  550. end do pole_loop4
  551. if (grid1_add /=0) then
  552. grid1_area(grid1_add) = grid1_area(grid1_add) + weights(1)
  553. grid1_centroid_lat(grid1_add) =
  554. & grid1_centroid_lat(grid1_add) + weights(2)
  555. grid1_centroid_lon(grid1_add) =
  556. & grid1_centroid_lon(grid1_add) + weights(3)
  557. endif
  558. if (grid2_add /=0) then
  559. grid2_area(grid2_add) = grid2_area(grid2_add) +
  560. & weights(num_wts+1)
  561. grid2_centroid_lat(grid2_add) =
  562. & grid2_centroid_lat(grid2_add) + weights(num_wts+2)
  563. grid2_centroid_lon(grid2_add) =
  564. & grid2_centroid_lon(grid2_add) + weights(num_wts+3)
  565. endif
  566. if (grid1_add /= 0 .and. grid2_add /=0) then
  567. call store_link_cnsrv(grid1_add, grid2_add, weights)
  568. grid1_frac(grid1_add) = grid1_frac(grid1_add) +
  569. & weights(1)
  570. grid2_frac(grid2_add) = grid2_frac(grid2_add) +
  571. & weights(num_wts+1)
  572. endif
  573. !-----------------------------------------------------------------------
  574. !
  575. ! finish centroid computation
  576. !
  577. !-----------------------------------------------------------------------
  578. where (grid1_area /= zero)
  579. grid1_centroid_lat = grid1_centroid_lat/grid1_area
  580. grid1_centroid_lon = grid1_centroid_lon/grid1_area
  581. end where
  582. where (grid2_area /= zero)
  583. grid2_centroid_lat = grid2_centroid_lat/grid2_area
  584. grid2_centroid_lon = grid2_centroid_lon/grid2_area
  585. end where
  586. !-----------------------------------------------------------------------
  587. !
  588. ! include centroids in weights and normalize using destination
  589. ! area if requested
  590. !
  591. !-----------------------------------------------------------------------
  592. do n=1,num_links_map1
  593. grid1_add = grid1_add_map1(n)
  594. grid2_add = grid2_add_map1(n)
  595. do nwgt=1,num_wts
  596. weights( nwgt) = wts_map1(nwgt,n)
  597. if (num_maps > 1) then
  598. weights(num_wts+nwgt) = wts_map2(nwgt,n)
  599. endif
  600. end do
  601. select case(norm_opt)
  602. case (norm_opt_dstarea)
  603. if (grid2_area(grid2_add) /= zero) then
  604. if (luse_grid2_area) then
  605. norm_factor = one/grid2_area_in(grid2_add)
  606. else
  607. norm_factor = one/grid2_area(grid2_add)
  608. endif
  609. else
  610. norm_factor = zero
  611. endif
  612. case (norm_opt_frcarea)
  613. if (grid2_frac(grid2_add) /= zero) then
  614. if (luse_grid2_area) then
  615. norm_factor = grid2_area(grid2_add)/
  616. & (grid2_frac(grid2_add)*
  617. & grid2_area_in(grid2_add))
  618. else
  619. norm_factor = one/grid2_frac(grid2_add)
  620. endif
  621. else
  622. norm_factor = zero
  623. endif
  624. case (norm_opt_none)
  625. norm_factor = one
  626. end select
  627. wts_map1(1,n) = weights(1)*norm_factor
  628. wts_map1(2,n) = (weights(2) - weights(1)*
  629. & grid1_centroid_lat(grid1_add))*
  630. & norm_factor
  631. wts_map1(3,n) = (weights(3) - weights(1)*
  632. & grid1_centroid_lon(grid1_add))*
  633. & norm_factor
  634. if (num_maps > 1) then
  635. select case(norm_opt)
  636. case (norm_opt_dstarea)
  637. if (grid1_area(grid1_add) /= zero) then
  638. if (luse_grid1_area) then
  639. norm_factor = one/grid1_area_in(grid1_add)
  640. else
  641. norm_factor = one/grid1_area(grid1_add)
  642. endif
  643. else
  644. norm_factor = zero
  645. endif
  646. case (norm_opt_frcarea)
  647. if (grid1_frac(grid1_add) /= zero) then
  648. if (luse_grid1_area) then
  649. norm_factor = grid1_area(grid1_add)/
  650. & (grid1_frac(grid1_add)*
  651. & grid1_area_in(grid1_add))
  652. else
  653. norm_factor = one/grid1_frac(grid1_add)
  654. endif
  655. else
  656. norm_factor = zero
  657. endif
  658. case (norm_opt_none)
  659. norm_factor = one
  660. end select
  661. wts_map2(1,n) = weights(num_wts+1)*norm_factor
  662. wts_map2(2,n) = (weights(num_wts+2) - weights(num_wts+1)*
  663. & grid2_centroid_lat(grid2_add))*
  664. & norm_factor
  665. wts_map2(3,n) = (weights(num_wts+3) - weights(num_wts+1)*
  666. & grid2_centroid_lon(grid2_add))*
  667. & norm_factor
  668. endif
  669. end do
  670. print *, 'Total number of links = ',num_links_map1
  671. where (grid1_area /= zero) grid1_frac = grid1_frac/grid1_area
  672. where (grid2_area /= zero) grid2_frac = grid2_frac/grid2_area
  673. !-----------------------------------------------------------------------
  674. !
  675. ! perform some error checking on final weights
  676. !
  677. !-----------------------------------------------------------------------
  678. grid2_centroid_lat = zero
  679. grid2_centroid_lon = zero
  680. do n=1,grid1_size
  681. if (grid1_area(n) < -.01) then
  682. print *,'Grid 1 area error: ',n,grid1_area(n)
  683. endif
  684. if (grid1_centroid_lat(n) < -pih-.01 .or.
  685. & grid1_centroid_lat(n) > pih+.01) then
  686. print *,'Grid 1 centroid lat error: ',n,grid1_centroid_lat(n)
  687. endif
  688. grid1_centroid_lat(n) = zero
  689. grid1_centroid_lon(n) = zero
  690. end do
  691. do n=1,grid2_size
  692. if (grid2_area(n) < -.01) then
  693. print *,'Grid 2 area error: ',n,grid2_area(n)
  694. endif
  695. if (grid2_centroid_lat(n) < -pih-.01 .or.
  696. & grid2_centroid_lat(n) > pih+.01) then
  697. print *,'Grid 2 centroid lat error: ',n,grid2_centroid_lat(n)
  698. endif
  699. grid2_centroid_lat(n) = zero
  700. grid2_centroid_lon(n) = zero
  701. end do
  702. do n=1,num_links_map1
  703. grid1_add = grid1_add_map1(n)
  704. grid2_add = grid2_add_map1(n)
  705. if (wts_map1(1,n) < -.01) then
  706. print *,'Map 1 weight < 0 ',grid1_add,grid2_add,wts_map1(1,n)
  707. endif
  708. if (norm_opt /= norm_opt_none .and. wts_map1(1,n) > 1.01) then
  709. print *,'Map 1 weight > 1 ',grid1_add,grid2_add,wts_map1(1,n)
  710. endif
  711. grid2_centroid_lat(grid2_add) =
  712. & grid2_centroid_lat(grid2_add) + wts_map1(1,n)
  713. if (num_maps > 1) then
  714. if (wts_map2(1,n) < -.01) then
  715. print *,'Map 2 weight < 0 ',grid1_add,grid2_add,
  716. & wts_map2(1,n)
  717. endif
  718. if (norm_opt /= norm_opt_none .and. wts_map2(1,n) > 1.01) then
  719. print *,'Map 2 weight < 0 ',grid1_add,grid2_add,
  720. & wts_map2(1,n)
  721. endif
  722. grid1_centroid_lat(grid1_add) =
  723. & grid1_centroid_lat(grid1_add) + wts_map2(1,n)
  724. endif
  725. end do
  726. do n=1,grid2_size
  727. select case(norm_opt)
  728. case (norm_opt_dstarea)
  729. norm_factor = grid2_frac(grid2_add)
  730. case (norm_opt_frcarea)
  731. norm_factor = one
  732. case (norm_opt_none)
  733. if (luse_grid2_area) then
  734. norm_factor = grid2_area_in(grid2_add)
  735. else
  736. norm_factor = grid2_area(grid2_add)
  737. endif
  738. end select
  739. if (abs(grid2_centroid_lat(grid2_add)-norm_factor) > .01) then
  740. print *,'Error: sum of wts for map1 ',grid2_add,
  741. & grid2_centroid_lat(grid2_add),norm_factor
  742. endif
  743. end do
  744. if (num_maps > 1) then
  745. do n=1,grid1_size
  746. select case(norm_opt)
  747. case (norm_opt_dstarea)
  748. norm_factor = grid1_frac(grid1_add)
  749. case (norm_opt_frcarea)
  750. norm_factor = one
  751. case (norm_opt_none)
  752. if (luse_grid1_area) then
  753. norm_factor = grid1_area_in(grid1_add)
  754. else
  755. norm_factor = grid1_area(grid1_add)
  756. endif
  757. end select
  758. if (abs(grid1_centroid_lat(grid1_add)-norm_factor) > .01) then
  759. print *,'Error: sum of wts for map2 ',grid1_add,
  760. & grid1_centroid_lat(grid1_add),norm_factor
  761. endif
  762. end do
  763. endif
  764. !-----------------------------------------------------------------------
  765. end subroutine remap_conserv
  766. !***********************************************************************
  767. subroutine intersection(location,intrsct_lat,intrsct_lon,lcoinc,
  768. & beglat, beglon, endlat, endlon, begseg,
  769. & lbegin, lrevers)
  770. !-----------------------------------------------------------------------
  771. !
  772. ! this routine finds the next intersection of a destination grid
  773. ! line with the line segment given by beglon, endlon, etc.
  774. ! a coincidence flag is returned if the segment is entirely
  775. ! coincident with an ocean grid line. the cells in which to search
  776. ! for an intersection must have already been restricted in the
  777. ! calling routine.
  778. !
  779. !-----------------------------------------------------------------------
  780. !-----------------------------------------------------------------------
  781. !
  782. ! intent(in):
  783. !
  784. !-----------------------------------------------------------------------
  785. logical (kind=log_kind), intent(in) ::
  786. & lbegin, ! flag for first integration along this segment
  787. & lrevers ! flag whether segment integrated in reverse
  788. real (kind=dbl_kind), intent(in) ::
  789. & beglat, beglon, ! beginning lat/lon endpoints for segment
  790. & endlat, endlon ! ending lat/lon endpoints for segment
  791. real (kind=dbl_kind), dimension(2), intent(inout) ::
  792. & begseg ! begin lat/lon of full segment
  793. !-----------------------------------------------------------------------
  794. !
  795. ! intent(out):
  796. !
  797. !-----------------------------------------------------------------------
  798. integer (kind=int_kind), intent(out) ::
  799. & location ! address in destination array containing this
  800. ! segment
  801. logical (kind=log_kind), intent(out) ::
  802. & lcoinc ! flag segments which are entirely coincident
  803. ! with a grid line
  804. real (kind=dbl_kind), intent(out) ::
  805. & intrsct_lat, intrsct_lon ! lat/lon coords of next intersect.
  806. !-----------------------------------------------------------------------
  807. !
  808. ! local variables
  809. !
  810. !-----------------------------------------------------------------------
  811. integer (kind=int_kind) :: n, next_n, cell, srch_corners, pole_loc
  812. integer (kind=int_kind), save ::
  813. & last_loc ! save location when crossing threshold
  814. logical (kind=log_kind) ::
  815. & loutside ! flags points outside grid
  816. logical (kind=log_kind), save ::
  817. & lthresh = .false. ! flags segments crossing threshold bndy
  818. real (kind=dbl_kind) ::
  819. & lon1, lon2, ! local longitude variables for segment
  820. & lat1, lat2, ! local latitude variables for segment
  821. & grdlon1, grdlon2, ! local longitude variables for grid cell
  822. & grdlat1, grdlat2, ! local latitude variables for grid cell
  823. & vec1_lat, vec1_lon, ! vectors and cross products used
  824. & vec2_lat, vec2_lon, ! during grid search
  825. & cross_product,
  826. & eps, offset, ! small offset away from intersect
  827. & s1, s2, determ, ! variables used for linear solve to
  828. & mat1, mat2, mat3, mat4, rhs1, rhs2 ! find intersection
  829. real (kind=dbl_kind), save ::
  830. & intrsct_lat_off, intrsct_lon_off ! lat/lon coords offset
  831. ! for next search
  832. !-----------------------------------------------------------------------
  833. !
  834. ! initialize defaults, flags, etc.
  835. !
  836. !-----------------------------------------------------------------------
  837. location = 0
  838. lcoinc = .false.
  839. intrsct_lat = endlat
  840. intrsct_lon = endlon
  841. if (num_srch_cells == 0) return
  842. if (beglat > north_thresh .or. beglat < south_thresh) then
  843. if (lthresh) location = last_loc
  844. call pole_intersection(location,
  845. & intrsct_lat,intrsct_lon,lcoinc,lthresh,
  846. & beglat, beglon, endlat, endlon, begseg, lrevers)
  847. if (lthresh) then
  848. last_loc = location
  849. intrsct_lat_off = intrsct_lat
  850. intrsct_lon_off = intrsct_lon
  851. endif
  852. return
  853. endif
  854. loutside = .false.
  855. if (lbegin) then
  856. lat1 = beglat
  857. lon1 = beglon
  858. else
  859. lat1 = intrsct_lat_off
  860. lon1 = intrsct_lon_off
  861. endif
  862. lat2 = endlat
  863. lon2 = endlon
  864. if ((lon2-lon1) > three*pih) then
  865. lon2 = lon2 - pi2
  866. else if ((lon2-lon1) < -three*pih) then
  867. lon2 = lon2 + pi2
  868. endif
  869. s1 = zero
  870. !-----------------------------------------------------------------------
  871. !
  872. ! search for location of this segment in ocean grid using cross
  873. ! product method to determine whether a point is enclosed by a cell
  874. !
  875. !-----------------------------------------------------------------------
  876. call timer_start(12)
  877. srch_corners = size(srch_corner_lat,DIM=1)
  878. srch_loop: do
  879. !***
  880. !*** if last segment crossed threshold, use that location
  881. !***
  882. if (lthresh) then
  883. do cell=1,num_srch_cells
  884. if (srch_add(cell) == last_loc) then
  885. location = last_loc
  886. eps = tiny
  887. exit srch_loop
  888. endif
  889. end do
  890. endif
  891. !***
  892. !*** otherwise normal search algorithm
  893. !***
  894. cell_loop: do cell=1,num_srch_cells
  895. corner_loop: do n=1,srch_corners
  896. next_n = MOD(n,srch_corners) + 1
  897. !***
  898. !*** here we take the cross product of the vector making
  899. !*** up each cell side with the vector formed by the vertex
  900. !*** and search point. if all the cross products are
  901. !*** positive, the point is contained in the cell.
  902. !***
  903. vec1_lat = srch_corner_lat(next_n,cell) -
  904. & srch_corner_lat(n ,cell)
  905. vec1_lon = srch_corner_lon(next_n,cell) -
  906. & srch_corner_lon(n ,cell)
  907. vec2_lat = lat1 - srch_corner_lat(n,cell)
  908. vec2_lon = lon1 - srch_corner_lon(n,cell)
  909. !***
  910. !*** if endpoint coincident with vertex, offset
  911. !*** the endpoint
  912. !***
  913. if (vec2_lat == 0 .and. vec2_lon == 0) then
  914. lat1 = lat1 + 1.d-10*(lat2-lat1)
  915. lon1 = lon1 + 1.d-10*(lon2-lon1)
  916. vec2_lat = lat1 - srch_corner_lat(n,cell)
  917. vec2_lon = lon1 - srch_corner_lon(n,cell)
  918. endif
  919. !***
  920. !*** check for 0,2pi crossings
  921. !***
  922. if (vec1_lon > pi) then
  923. vec1_lon = vec1_lon - pi2
  924. else if (vec1_lon < -pi) then
  925. vec1_lon = vec1_lon + pi2
  926. endif
  927. if (vec2_lon > pi) then
  928. vec2_lon = vec2_lon - pi2
  929. else if (vec2_lon < -pi) then
  930. vec2_lon = vec2_lon + pi2
  931. endif
  932. cross_product = vec1_lon*vec2_lat - vec2_lon*vec1_lat
  933. !***
  934. !*** if the cross product for a side is zero, the point
  935. !*** lies exactly on the side or the side is degenerate
  936. !*** (zero length). if degenerate, set the cross
  937. !*** product to a positive number. otherwise perform
  938. !*** another cross product between the side and the
  939. !*** segment itself.
  940. !*** if this cross product is also zero, the line is
  941. !*** coincident with the cell boundary - perform the
  942. !*** dot product and only choose the cell if the dot
  943. !*** product is positive (parallel vs anti-parallel).
  944. !***
  945. if (cross_product == zero) then
  946. if (vec1_lat /= zero .or. vec1_lon /= zero) then
  947. vec2_lat = lat2 - lat1
  948. vec2_lon = lon2 - lon1
  949. if (vec2_lon > pi) then
  950. vec2_lon = vec2_lon - pi2
  951. else if (vec2_lon < -pi) then
  952. vec2_lon = vec2_lon + pi2
  953. endif
  954. cross_product = vec1_lon*vec2_lat - vec2_lon*vec1_lat
  955. else
  956. cross_product = one
  957. endif
  958. if (cross_product == zero) then
  959. lcoinc = .true.
  960. cross_product = vec1_lon*vec2_lon + vec1_lat*vec2_lat
  961. if (lrevers) cross_product = -cross_product
  962. endif
  963. endif
  964. !***
  965. !*** if cross product is less than zero, this cell
  966. !*** doesn't work
  967. !***
  968. if (cross_product < zero) exit corner_loop
  969. end do corner_loop
  970. !***
  971. !*** if cross products all positive, we found the location
  972. !***
  973. if (n > srch_corners) then
  974. location = srch_add(cell)
  975. !***
  976. !*** if the beginning of this segment was outside the
  977. !*** grid, invert the segment so the intersection found
  978. !*** will be the first intersection with the grid
  979. !***
  980. if (loutside) then
  981. lat2 = beglat
  982. lon2 = beglon
  983. location = 0
  984. eps = -tiny
  985. else
  986. eps = tiny
  987. endif
  988. exit srch_loop
  989. endif
  990. !***
  991. !*** otherwise move on to next cell
  992. !***
  993. end do cell_loop
  994. !***
  995. !*** if still no cell found, the point lies outside the grid.
  996. !*** take some baby steps along the segment to see if any
  997. !*** part of the segment lies inside the grid.
  998. !***
  999. loutside = .true.
  1000. s1 = s1 + 0.001_dbl_kind
  1001. lat1 = beglat + s1*(endlat - beglat)
  1002. lon1 = beglon + s1*(lon2 - beglon)
  1003. !***
  1004. !*** reached the end of the segment and still outside the grid
  1005. !*** return no intersection
  1006. !***
  1007. if (s1 >= one) return
  1008. end do srch_loop
  1009. call timer_stop(12)
  1010. !-----------------------------------------------------------------------
  1011. !
  1012. ! now that a cell is found, search for the next intersection.
  1013. ! loop over sides of the cell to find intersection with side
  1014. ! must check all sides for coincidences or intersections
  1015. !
  1016. !-----------------------------------------------------------------------
  1017. call timer_start(13)
  1018. intrsct_loop: do n=1,srch_corners
  1019. next_n = mod(n,srch_corners) + 1
  1020. grdlon1 = srch_corner_lon(n ,cell)
  1021. grdlon2 = srch_corner_lon(next_n,cell)
  1022. grdlat1 = srch_corner_lat(n ,cell)
  1023. grdlat2 = srch_corner_lat(next_n,cell)
  1024. !***
  1025. !*** set up linear system to solve for intersection
  1026. !***
  1027. mat1 = lat2 - lat1
  1028. mat2 = grdlat1 - grdlat2
  1029. mat3 = lon2 - lon1
  1030. mat4 = grdlon1 - grdlon2
  1031. rhs1 = grdlat1 - lat1
  1032. rhs2 = grdlon1 - lon1
  1033. if (mat3 > pi) then
  1034. mat3 = mat3 - pi2
  1035. else if (mat3 < -pi) then
  1036. mat3 = mat3 + pi2
  1037. endif
  1038. if (mat4 > pi) then
  1039. mat4 = mat4 - pi2
  1040. else if (mat4 < -pi) then
  1041. mat4 = mat4 + pi2
  1042. endif
  1043. if (rhs2 > pi) then
  1044. rhs2 = rhs2 - pi2
  1045. else if (rhs2 < -pi) then
  1046. rhs2 = rhs2 + pi2
  1047. endif
  1048. determ = mat1*mat4 - mat2*mat3
  1049. !***
  1050. !*** if the determinant is zero, the segments are either
  1051. !*** parallel or coincident. coincidences were detected
  1052. !*** above so do nothing.
  1053. !*** if the determinant is non-zero, solve for the linear
  1054. !*** parameters s for the intersection point on each line
  1055. !*** segment.
  1056. !*** if 0<s1,s2<1 then the segment intersects with this side.
  1057. !*** return the point of intersection (adding a small
  1058. !*** number so the intersection is off the grid line).
  1059. !***
  1060. if (abs(determ) > 1.e-30) then
  1061. s1 = (rhs1*mat4 - mat2*rhs2)/determ
  1062. s2 = (mat1*rhs2 - rhs1*mat3)/determ
  1063. if (s2 >= zero .and. s2 <= one .and.
  1064. & s1 > zero. and. s1 <= one) then
  1065. !***
  1066. !*** recompute intersection based on full segment
  1067. !*** so intersections are consistent for both sweeps
  1068. !***
  1069. if (.not. loutside) then
  1070. mat1 = lat2 - begseg(1)
  1071. mat3 = lon2 - begseg(2)
  1072. rhs1 = grdlat1 - begseg(1)
  1073. rhs2 = grdlon1 - begseg(2)
  1074. else
  1075. mat1 = begseg(1) - endlat
  1076. mat3 = begseg(2) - endlon
  1077. rhs1 = grdlat1 - endlat
  1078. rhs2 = grdlon1 - endlon
  1079. endif
  1080. if (mat3 > pi) then
  1081. mat3 = mat3 - pi2
  1082. else if (mat3 < -pi) then
  1083. mat3 = mat3 + pi2
  1084. endif
  1085. if (rhs2 > pi) then
  1086. rhs2 = rhs2 - pi2
  1087. else if (rhs2 < -pi) then
  1088. rhs2 = rhs2 + pi2
  1089. endif
  1090. determ = mat1*mat4 - mat2*mat3
  1091. !***
  1092. !*** sometimes due to roundoff, the previous
  1093. !*** determinant is non-zero, but the lines
  1094. !*** are actually coincident. if this is the
  1095. !*** case, skip the rest.
  1096. !***
  1097. if (determ /= zero) then
  1098. s1 = (rhs1*mat4 - mat2*rhs2)/determ
  1099. s2 = (mat1*rhs2 - rhs1*mat3)/determ
  1100. offset = s1 + eps/determ
  1101. if (offset > one) offset = one
  1102. if (.not. loutside) then
  1103. intrsct_lat = begseg(1) + mat1*s1
  1104. intrsct_lon = begseg(2) + mat3*s1
  1105. intrsct_lat_off = begseg(1) + mat1*offset
  1106. intrsct_lon_off = begseg(2) + mat3*offset
  1107. else
  1108. intrsct_lat = endlat + mat1*s1
  1109. intrsct_lon = endlon + mat3*s1
  1110. intrsct_lat_off = endlat + mat1*offset
  1111. intrsct_lon_off = endlon + mat3*offset
  1112. endif
  1113. exit intrsct_loop
  1114. endif
  1115. endif
  1116. endif
  1117. !***
  1118. !*** no intersection this side, move on to next side
  1119. !***
  1120. end do intrsct_loop
  1121. call timer_stop(13)
  1122. !-----------------------------------------------------------------------
  1123. !
  1124. ! if the segment crosses a pole threshold, reset the intersection
  1125. ! to be the threshold latitude. only check if this was not a
  1126. ! threshold segment since sometimes coordinate transform can end
  1127. ! up on other side of threshold again.
  1128. !
  1129. !-----------------------------------------------------------------------
  1130. if (lthresh) then
  1131. if (intrsct_lat < north_thresh .or. intrsct_lat > south_thresh)
  1132. & lthresh = .false.
  1133. else if (lat1 > zero .and. intrsct_lat > north_thresh) then
  1134. intrsct_lat = north_thresh + tiny
  1135. intrsct_lat_off = north_thresh + eps*mat1
  1136. s1 = (intrsct_lat - begseg(1))/mat1
  1137. intrsct_lon = begseg(2) + s1*mat3
  1138. intrsct_lon_off = begseg(2) + (s1+eps)*mat3
  1139. last_loc = location
  1140. lthresh = .true.
  1141. else if (lat1 < zero .and. intrsct_lat < south_thresh) then
  1142. intrsct_lat = south_thresh - tiny
  1143. intrsct_lat_off = south_thresh + eps*mat1
  1144. s1 = (intrsct_lat - begseg(1))/mat1
  1145. intrsct_lon = begseg(2) + s1*mat3
  1146. intrsct_lon_off = begseg(2) + (s1+eps)*mat3
  1147. last_loc = location
  1148. lthresh = .true.
  1149. endif
  1150. !-----------------------------------------------------------------------
  1151. end subroutine intersection
  1152. !***********************************************************************
  1153. subroutine pole_intersection(location,
  1154. & intrsct_lat,intrsct_lon,lcoinc,lthresh,
  1155. & beglat, beglon, endlat, endlon, begseg, lrevers)
  1156. !-----------------------------------------------------------------------
  1157. !
  1158. ! this routine is identical to the intersection routine except
  1159. ! that a coordinate transformation (using a Lambert azimuthal
  1160. ! equivalent projection) is performed to treat polar cells more
  1161. ! accurately.
  1162. !
  1163. !-----------------------------------------------------------------------
  1164. !-----------------------------------------------------------------------
  1165. !
  1166. ! intent(in):
  1167. !
  1168. !-----------------------------------------------------------------------
  1169. real (kind=dbl_kind), intent(in) ::
  1170. & beglat, beglon, ! beginning lat/lon endpoints for segment
  1171. & endlat, endlon ! ending lat/lon endpoints for segment
  1172. real (kind=dbl_kind), dimension(2), intent(inout) ::
  1173. & begseg ! begin lat/lon of full segment
  1174. logical (kind=log_kind), intent(in) ::
  1175. & lrevers ! flag true if segment integrated in reverse
  1176. !-----------------------------------------------------------------------
  1177. !
  1178. ! intent(out):
  1179. !
  1180. !-----------------------------------------------------------------------
  1181. integer (kind=int_kind), intent(inout) ::
  1182. & location ! address in destination array containing this
  1183. ! segment -- also may contain last location on
  1184. ! entry
  1185. logical (kind=log_kind), intent(out) ::
  1186. & lcoinc ! flag segment coincident with grid line
  1187. logical (kind=log_kind), intent(inout) ::
  1188. & lthresh ! flag segment crossing threshold boundary
  1189. real (kind=dbl_kind), intent(out) ::
  1190. & intrsct_lat, intrsct_lon ! lat/lon coords of next intersect.
  1191. !-----------------------------------------------------------------------
  1192. !
  1193. ! local variables
  1194. !
  1195. !-----------------------------------------------------------------------
  1196. integer (kind=int_kind) :: n, next_n, cell, srch_corners, pole_loc
  1197. logical (kind=log_kind) :: loutside ! flags points outside grid
  1198. real (kind=dbl_kind) :: pi4, rns, ! north/south conversion
  1199. & x1, x2, ! local x variables for segment
  1200. & y1, y2, ! local y variables for segment
  1201. & begx, begy, ! beginning x,y variables for segment
  1202. & endx, endy, ! beginning x,y variables for segment
  1203. & begsegx, begsegy, ! beginning x,y variables for segment
  1204. & grdx1, grdx2, ! local x variables for grid cell
  1205. & grdy1, grdy2, ! local y variables for grid cell
  1206. & vec1_y, vec1_x, ! vectors and cross products used
  1207. & vec2_y, vec2_x, ! during grid search
  1208. & cross_product, eps, ! eps=small offset away from intersect
  1209. & s1, s2, determ, ! variables used for linear solve to
  1210. & mat1, mat2, mat3, mat4, rhs1, rhs2 ! find intersection
  1211. real (kind=dbl_kind), dimension(:,:), allocatable ::
  1212. & srch_corner_x, ! x of each corner of srch cells
  1213. & srch_corner_y ! y of each corner of srch cells
  1214. !***
  1215. !*** save last intersection to avoid roundoff during coord
  1216. !*** transformation
  1217. !***
  1218. logical (kind=log_kind), save :: luse_last = .false.
  1219. real (kind=dbl_kind), save ::
  1220. & intrsct_x, intrsct_y ! x,y for intersection
  1221. !***
  1222. !*** variables necessary if segment manages to hit pole
  1223. !***
  1224. integer (kind=int_kind), save ::
  1225. & avoid_pole_count = 0 ! count attempts to avoid pole
  1226. real (kind=dbl_kind), save ::
  1227. & avoid_pole_offset = tiny ! endpoint offset to avoid pole
  1228. !-----------------------------------------------------------------------
  1229. !
  1230. ! initialize defaults, flags, etc.
  1231. !
  1232. !-----------------------------------------------------------------------
  1233. if (.not. lthresh) location = 0
  1234. lcoinc = .false.
  1235. intrsct_lat = endlat
  1236. intrsct_lon = endlon
  1237. loutside = .false.
  1238. s1 = zero
  1239. !-----------------------------------------------------------------------
  1240. !
  1241. ! convert coordinates
  1242. !
  1243. !-----------------------------------------------------------------------
  1244. allocate(srch_corner_x(size(srch_corner_lat,DIM=1),
  1245. & size(srch_corner_lat,DIM=2)),
  1246. & srch_corner_y(size(srch_corner_lat,DIM=1),
  1247. & size(srch_corner_lat,DIM=2)))
  1248. if (beglat > zero) then
  1249. pi4 = quart*pi
  1250. rns = one
  1251. else
  1252. pi4 = -quart*pi
  1253. rns = -one
  1254. endif
  1255. if (luse_last) then
  1256. x1 = intrsct_x
  1257. y1 = intrsct_y
  1258. else
  1259. x1 = rns*two*sin(pi4 - half*beglat)*cos(beglon)
  1260. y1 = two*sin(pi4 - half*beglat)*sin(beglon)
  1261. luse_last = .true.
  1262. endif
  1263. x2 = rns*two*sin(pi4 - half*endlat)*cos(endlon)
  1264. y2 = two*sin(pi4 - half*endlat)*sin(endlon)
  1265. srch_corner_x = rns*two*sin(pi4 - half*srch_corner_lat)*
  1266. & cos(srch_corner_lon)
  1267. srch_corner_y = two*sin(pi4 - half*srch_corner_lat)*
  1268. & sin(srch_corner_lon)
  1269. begx = x1
  1270. begy = y1
  1271. endx = x2
  1272. endy = y2
  1273. begsegx = rns*two*sin(pi4 - half*begseg(1))*cos(begseg(2))
  1274. begsegy = two*sin(pi4 - half*begseg(1))*sin(begseg(2))
  1275. intrsct_x = endx
  1276. intrsct_y = endy
  1277. !-----------------------------------------------------------------------
  1278. !
  1279. ! search for location of this segment in ocean grid using cross
  1280. ! product method to determine whether a point is enclosed by a cell
  1281. !
  1282. !-----------------------------------------------------------------------
  1283. call timer_start(12)
  1284. srch_corners = size(srch_corner_lat,DIM=1)
  1285. srch_loop: do
  1286. !***
  1287. !*** if last segment crossed threshold, use that location
  1288. !***
  1289. if (lthresh) then
  1290. do cell=1,num_srch_cells
  1291. if (srch_add(cell) == location) then
  1292. eps = tiny
  1293. exit srch_loop
  1294. endif
  1295. end do
  1296. endif
  1297. !***
  1298. !*** otherwise normal search algorithm
  1299. !***
  1300. cell_loop: do cell=1,num_srch_cells
  1301. corner_loop: do n=1,srch_corners
  1302. next_n = MOD(n,srch_corners) + 1
  1303. !***
  1304. !*** here we take the cross product of the vector making
  1305. !*** up each cell side with the vector formed by the vertex
  1306. !*** and search point. if all the cross products are
  1307. !*** positive, the point is contained in the cell.
  1308. !***
  1309. vec1_x = srch_corner_x(next_n,cell) -
  1310. & srch_corner_x(n ,cell)
  1311. vec1_y = srch_corner_y(next_n,cell) -
  1312. & srch_corner_y(n ,cell)
  1313. vec2_x = x1 - srch_corner_x(n,cell)
  1314. vec2_y = y1 - srch_corner_y(n,cell)
  1315. !***
  1316. !*** if endpoint coincident with vertex, offset
  1317. !*** the endpoint
  1318. !***
  1319. if (vec2_x == 0 .and. vec2_y == 0) then
  1320. x1 = x1 + 1.d-10*(x2-x1)
  1321. y1 = y1 + 1.d-10*(y2-y1)
  1322. vec2_x = x1 - srch_corner_x(n,cell)
  1323. vec2_y = y1 - srch_corner_y(n,cell)
  1324. endif
  1325. cross_product = vec1_x*vec2_y - vec2_x*vec1_y
  1326. !***
  1327. !*** if the cross product for a side is zero, the point
  1328. !*** lies exactly on the side or the length of a side
  1329. !*** is zero. if the length is zero set det > 0.
  1330. !*** otherwise, perform another cross
  1331. !*** product between the side and the segment itself.
  1332. !*** if this cross product is also zero, the line is
  1333. !*** coincident with the cell boundary - perform the
  1334. !*** dot product and only choose the cell if the dot
  1335. !*** product is positive (parallel vs anti-parallel).
  1336. !***
  1337. if (cross_product == zero) then
  1338. if (vec1_x /= zero .or. vec1_y /= 0) then
  1339. vec2_x = x2 - x1
  1340. vec2_y = y2 - y1
  1341. cross_product = vec1_x*vec2_y - vec2_x*vec1_y
  1342. else
  1343. cross_product = one
  1344. endif
  1345. if (cross_product == zero) then
  1346. lcoinc = .true.
  1347. cross_product = vec1_x*vec2_x + vec1_y*vec2_y
  1348. if (lrevers) cross_product = -cross_product
  1349. endif
  1350. endif
  1351. !***
  1352. !*** if cross product is less than zero, this cell
  1353. !*** doesn't work
  1354. !***
  1355. if (cross_product < zero) exit corner_loop
  1356. end do corner_loop
  1357. !***
  1358. !*** if cross products all positive, we found the location
  1359. !***
  1360. if (n > srch_corners) then
  1361. location = srch_add(cell)
  1362. !***
  1363. !*** if the beginning of this segment was outside the
  1364. !*** grid, invert the segment so the intersection found
  1365. !*** will be the first intersection with the grid
  1366. !***
  1367. if (loutside) then
  1368. x2 = begx
  1369. y2 = begy
  1370. location = 0
  1371. eps = -tiny
  1372. else
  1373. eps = tiny
  1374. endif
  1375. exit srch_loop
  1376. endif
  1377. !***
  1378. !*** otherwise move on to next cell
  1379. !***
  1380. end do cell_loop
  1381. !***
  1382. !*** if no cell found, the point lies outside the grid.
  1383. !*** take some baby steps along the segment to see if any
  1384. !*** part of the segment lies inside the grid.
  1385. !***
  1386. loutside = .true.
  1387. s1 = s1 + 0.001_dbl_kind
  1388. x1 = begx + s1*(x2 - begx)
  1389. y1 = begy + s1*(y2 - begy)
  1390. !***
  1391. !*** reached the end of the segment and still outside the grid
  1392. !*** return no intersection
  1393. !***
  1394. if (s1 >= one) then
  1395. deallocate(srch_corner_x, srch_corner_y)
  1396. luse_last = .false.
  1397. return
  1398. endif
  1399. end do srch_loop
  1400. call timer_stop(12)
  1401. !-----------------------------------------------------------------------
  1402. !
  1403. ! now that a cell is found, search for the next intersection.
  1404. ! loop over sides of the cell to find intersection with side
  1405. ! must check all sides for coincidences or intersections
  1406. !
  1407. !-----------------------------------------------------------------------
  1408. call timer_start(13)
  1409. intrsct_loop: do n=1,srch_corners
  1410. next_n = mod(n,srch_corners) + 1
  1411. grdy1 = srch_corner_y(n ,cell)
  1412. grdy2 = srch_corner_y(next_n,cell)
  1413. grdx1 = srch_corner_x(n ,cell)
  1414. grdx2 = srch_corner_x(next_n,cell)
  1415. !***
  1416. !*** set up linear system to solve for intersection
  1417. !***
  1418. mat1 = x2 - x1
  1419. mat2 = grdx1 - grdx2
  1420. mat3 = y2 - y1
  1421. mat4 = grdy1 - grdy2
  1422. rhs1 = grdx1 - x1
  1423. rhs2 = grdy1 - y1
  1424. determ = mat1*mat4 - mat2*mat3
  1425. !***
  1426. !*** if the determinant is zero, the segments are either
  1427. !*** parallel or coincident or one segment has zero length.
  1428. !*** coincidences were detected above so do nothing.
  1429. !*** if the determinant is non-zero, solve for the linear
  1430. !*** parameters s for the intersection point on each line
  1431. !*** segment.
  1432. !*** if 0<s1,s2<1 then the segment intersects with this side.
  1433. !*** return the point of intersection (adding a small
  1434. !*** number so the intersection is off the grid line).
  1435. !***
  1436. if (abs(determ) > 1.e-30) then
  1437. s1 = (rhs1*mat4 - mat2*rhs2)/determ
  1438. s2 = (mat1*rhs2 - rhs1*mat3)/determ
  1439. if (s2 >= zero .and. s2 <= one .and.
  1440. & s1 > zero. and. s1 <= one) then
  1441. !***
  1442. !*** recompute intersection using entire segment
  1443. !*** for consistency between sweeps
  1444. !***
  1445. if (.not. loutside) then
  1446. mat1 = x2 - begsegx
  1447. mat3 = y2 - begsegy
  1448. rhs1 = grdx1 - begsegx
  1449. rhs2 = grdy1 - begsegy
  1450. else
  1451. mat1 = x2 - endx
  1452. mat3 = y2 - endy
  1453. rhs1 = grdx1 - endx
  1454. rhs2 = grdy1 - endy
  1455. endif
  1456. determ = mat1*mat4 - mat2*mat3
  1457. !***
  1458. !*** sometimes due to roundoff, the previous
  1459. !*** determinant is non-zero, but the lines
  1460. !*** are actually coincident. if this is the
  1461. !*** case, skip the rest.
  1462. !***
  1463. if (determ /= zero) then
  1464. s1 = (rhs1*mat4 - mat2*rhs2)/determ
  1465. s2 = (mat1*rhs2 - rhs1*mat3)/determ
  1466. if (.not. loutside) then
  1467. intrsct_x = begsegx + s1*mat1
  1468. intrsct_y = begsegy + s1*mat3
  1469. else
  1470. intrsct_x = endx + s1*mat1
  1471. intrsct_y = endy + s1*mat3
  1472. endif
  1473. !***
  1474. !*** convert back to lat/lon coordinates
  1475. !***
  1476. intrsct_lon = rns*atan2(intrsct_y,intrsct_x)
  1477. if (intrsct_lon < zero)
  1478. & intrsct_lon = intrsct_lon + pi2
  1479. if (abs(intrsct_x) > 1.d-10) then
  1480. intrsct_lat = (pi4 -
  1481. & asin(rns*half*intrsct_x/cos(intrsct_lon)))*two
  1482. else if (abs(intrsct_y) > 1.d-10) then
  1483. intrsct_lat = (pi4 -
  1484. & asin(half*intrsct_y/sin(intrsct_lon)))*two
  1485. else
  1486. intrsct_lat = two*pi4
  1487. endif
  1488. !***
  1489. !*** add offset in transformed space for next pass.
  1490. !***
  1491. if (s1 - eps/determ < one) then
  1492. intrsct_x = intrsct_x - mat1*(eps/determ)
  1493. intrsct_y = intrsct_y - mat3*(eps/determ)
  1494. else
  1495. if (.not. loutside) then
  1496. intrsct_x = endx
  1497. intrsct_y = endy
  1498. intrsct_lat = endlat
  1499. intrsct_lon = endlon
  1500. else
  1501. intrsct_x = begsegx
  1502. intrsct_y = begsegy
  1503. intrsct_lat = begseg(1)
  1504. intrsct_lon = begseg(2)
  1505. endif
  1506. endif
  1507. exit intrsct_loop
  1508. endif
  1509. endif
  1510. endif
  1511. !***
  1512. !*** no intersection this side, move on to next side
  1513. !***
  1514. end do intrsct_loop
  1515. call timer_stop(13)
  1516. deallocate(srch_corner_x, srch_corner_y)
  1517. !-----------------------------------------------------------------------
  1518. !
  1519. ! if segment manages to cross over pole, shift the beginning
  1520. ! endpoint in order to avoid hitting pole directly
  1521. ! (it is ok for endpoint to be pole point)
  1522. !
  1523. !-----------------------------------------------------------------------
  1524. if (abs(intrsct_x) < 1.e-10 .and. abs(intrsct_y) < 1.e-10 .and.
  1525. & (endx /= zero .and. endy /=0)) then
  1526. if (avoid_pole_count > 2) then
  1527. avoid_pole_count = 0
  1528. avoid_pole_offset = 10.*avoid_pole_offset
  1529. endif
  1530. cross_product = begsegx*(endy-begsegy) - begsegy*(endx-begsegx)
  1531. intrsct_lat = begseg(1)
  1532. if (cross_product*intrsct_lat > zero) then
  1533. intrsct_lon = beglon + avoid_pole_offset
  1534. begseg(2) = begseg(2) + avoid_pole_offset
  1535. else
  1536. intrsct_lon = beglon - avoid_pole_offset
  1537. begseg(2) = begseg(2) - avoid_pole_offset
  1538. endif
  1539. avoid_pole_count = avoid_pole_count + 1
  1540. luse_last = .false.
  1541. else
  1542. avoid_pole_count = 0
  1543. avoid_pole_offset = tiny
  1544. endif
  1545. !-----------------------------------------------------------------------
  1546. !
  1547. ! if the segment crosses a pole threshold, reset the intersection
  1548. ! to be the threshold latitude and do not reuse x,y intersect
  1549. ! on next entry. only check if did not cross threshold last
  1550. ! time - sometimes the coordinate transformation can place a
  1551. ! segment on the other side of the threshold again
  1552. !
  1553. !-----------------------------------------------------------------------
  1554. if (lthresh) then
  1555. if (intrsct_lat > north_thresh .or. intrsct_lat < south_thresh)
  1556. & lthresh = .false.
  1557. else if (beglat > zero .and. intrsct_lat < north_thresh) then
  1558. mat4 = endlat - begseg(1)
  1559. mat3 = endlon - begseg(2)
  1560. if (mat3 > pi) mat3 = mat3 - pi2
  1561. if (mat3 < -pi) mat3 = mat3 + pi2
  1562. intrsct_lat = north_thresh - tiny
  1563. s1 = (north_thresh - begseg(1))/mat4
  1564. intrsct_lon = begseg(2) + s1*mat3
  1565. luse_last = .false.
  1566. lthresh = .true.
  1567. else if (beglat < zero .and. intrsct_lat > south_thresh) then
  1568. mat4 = endlat - begseg(1)
  1569. mat3 = endlon - begseg(2)
  1570. if (mat3 > pi) mat3 = mat3 - pi2
  1571. if (mat3 < -pi) mat3 = mat3 + pi2
  1572. intrsct_lat = south_thresh + tiny
  1573. s1 = (south_thresh - begseg(1))/mat4
  1574. intrsct_lon = begseg(2) + s1*mat3
  1575. luse_last = .false.
  1576. lthresh = .true.
  1577. endif
  1578. !***
  1579. !*** if reached end of segment, do not use x,y intersect
  1580. !*** on next entry
  1581. !***
  1582. if (intrsct_lat == endlat .and. intrsct_lon == endlon) then
  1583. luse_last = .false.
  1584. endif
  1585. !-----------------------------------------------------------------------
  1586. end subroutine pole_intersection
  1587. !***********************************************************************
  1588. subroutine line_integral(weights, num_wts,
  1589. & in_phi1, in_phi2, theta1, theta2,
  1590. & grid1_lat, grid1_lon, grid2_lat, grid2_lon)
  1591. !-----------------------------------------------------------------------
  1592. !
  1593. ! this routine computes the line integral of the flux function
  1594. ! that results in the interpolation weights. the line is defined
  1595. ! by the input lat/lon of the endpoints.
  1596. !
  1597. !-----------------------------------------------------------------------
  1598. !-----------------------------------------------------------------------
  1599. !
  1600. ! intent(in):
  1601. !
  1602. !-----------------------------------------------------------------------
  1603. integer (kind=int_kind), intent(in) ::
  1604. & num_wts ! number of weights to compute
  1605. real (kind=dbl_kind), intent(in) ::
  1606. & in_phi1, in_phi2, ! longitude endpoints for the segment
  1607. & theta1, theta2, ! latitude endpoints for the segment
  1608. & grid1_lat, grid1_lon, ! reference coordinates for each
  1609. & grid2_lat, grid2_lon ! grid (to ensure correct 0,2pi interv.
  1610. !-----------------------------------------------------------------------
  1611. !
  1612. ! intent(out):
  1613. !
  1614. !-----------------------------------------------------------------------
  1615. real (kind=dbl_kind), dimension(2*num_wts), intent(out) ::
  1616. & weights ! line integral contribution to weights
  1617. !-----------------------------------------------------------------------
  1618. !
  1619. ! local variables
  1620. !
  1621. !-----------------------------------------------------------------------
  1622. real (kind=dbl_kind) :: dphi, sinth1, sinth2, costh1, costh2, fac,
  1623. & phi1, phi2, phidiff1, phidiff2, sinint
  1624. real (kind=dbl_kind) :: f1, f2, fint
  1625. !-----------------------------------------------------------------------
  1626. !
  1627. ! weights for the general case based on a trapezoidal approx to
  1628. ! the integrals.
  1629. !
  1630. !-----------------------------------------------------------------------
  1631. sinth1 = SIN(theta1)
  1632. sinth2 = SIN(theta2)
  1633. costh1 = COS(theta1)
  1634. costh2 = COS(theta2)
  1635. dphi = in_phi1 - in_phi2
  1636. if (dphi > pi) then
  1637. dphi = dphi - pi2
  1638. else if (dphi < -pi) then
  1639. dphi = dphi + pi2
  1640. endif
  1641. dphi = half*dphi
  1642. !-----------------------------------------------------------------------
  1643. !
  1644. ! the first weight is the area overlap integral. the second and
  1645. ! fourth are second-order latitude gradient weights.
  1646. !
  1647. !-----------------------------------------------------------------------
  1648. weights( 1) = dphi*(sinth1 + sinth2)
  1649. weights(num_wts+1) = dphi*(sinth1 + sinth2)
  1650. weights( 2) = dphi*(costh1 + costh2 + (theta1*sinth1 +
  1651. & theta2*sinth2))
  1652. weights(num_wts+2) = dphi*(costh1 + costh2 + (theta1*sinth1 +
  1653. & theta2*sinth2))
  1654. !-----------------------------------------------------------------------
  1655. !
  1656. ! the third and fifth weights are for the second-order phi gradient
  1657. ! component. must be careful of longitude range.
  1658. !
  1659. !-----------------------------------------------------------------------
  1660. f1 = half*(costh1*sinth1 + theta1)
  1661. f2 = half*(costh2*sinth2 + theta2)
  1662. phi1 = in_phi1 - grid1_lon
  1663. if (phi1 > pi) then
  1664. phi1 = phi1 - pi2
  1665. else if (phi1 < -pi) then
  1666. phi1 = phi1 + pi2
  1667. endif
  1668. phi2 = in_phi2 - grid1_lon
  1669. if (phi2 > pi) then
  1670. phi2 = phi2 - pi2
  1671. else if (phi2 < -pi) then
  1672. phi2 = phi2 + pi2
  1673. endif
  1674. if ((phi2-phi1) < pi .and. (phi2-phi1) > -pi) then
  1675. weights(3) = dphi*(phi1*f1 + phi2*f2)
  1676. else
  1677. if (phi1 > zero) then
  1678. fac = pi
  1679. else
  1680. fac = -pi
  1681. endif
  1682. fint = f1 + (f2-f1)*(fac-phi1)/abs(dphi)
  1683. weights(3) = half*phi1*(phi1-fac)*f1 -
  1684. & half*phi2*(phi2+fac)*f2 +
  1685. & half*fac*(phi1+phi2)*fint
  1686. endif
  1687. phi1 = in_phi1 - grid2_lon
  1688. if (phi1 > pi) then
  1689. phi1 = phi1 - pi2
  1690. else if (phi1 < -pi) then
  1691. phi1 = phi1 + pi2
  1692. endif
  1693. phi2 = in_phi2 - grid2_lon
  1694. if (phi2 > pi) then
  1695. phi2 = phi2 - pi2
  1696. else if (phi2 < -pi) then
  1697. phi2 = phi2 + pi2
  1698. endif
  1699. if ((phi2-phi1) < pi .and. (phi2-phi1) > -pi) then
  1700. weights(num_wts+3) = dphi*(phi1*f1 + phi2*f2)
  1701. else
  1702. if (phi1 > zero) then
  1703. fac = pi
  1704. else
  1705. fac = -pi
  1706. endif
  1707. fint = f1 + (f2-f1)*(fac-phi1)/abs(dphi)
  1708. weights(num_wts+3) = half*phi1*(phi1-fac)*f1 -
  1709. & half*phi2*(phi2+fac)*f2 +
  1710. & half*fac*(phi1+phi2)*fint
  1711. endif
  1712. !-----------------------------------------------------------------------
  1713. end subroutine line_integral
  1714. !***********************************************************************
  1715. subroutine store_link_cnsrv(add1, add2, weights)
  1716. !-----------------------------------------------------------------------
  1717. !
  1718. ! this routine stores the address and weight for this link in
  1719. ! the appropriate address and weight arrays and resizes those
  1720. ! arrays if necessary.
  1721. !
  1722. !-----------------------------------------------------------------------
  1723. !-----------------------------------------------------------------------
  1724. !
  1725. ! input variables
  1726. !
  1727. !-----------------------------------------------------------------------
  1728. integer (kind=int_kind), intent(in) ::
  1729. & add1, ! address on grid1
  1730. & add2 ! address on grid2
  1731. real (kind=dbl_kind), dimension(:), intent(in) ::
  1732. & weights ! array of remapping weights for this link
  1733. !-----------------------------------------------------------------------
  1734. !
  1735. ! local variables
  1736. !
  1737. !-----------------------------------------------------------------------
  1738. integer (kind=int_kind) :: nlink, min_link, max_link ! link index
  1739. integer (kind=int_kind), dimension(:,:), allocatable, save ::
  1740. & link_add1, ! min,max link add to restrict search
  1741. & link_add2 ! min,max link add to restrict search
  1742. logical (kind=log_kind), save :: first_call = .true.
  1743. !-----------------------------------------------------------------------
  1744. !
  1745. ! if all weights are zero, do not bother storing the link
  1746. !
  1747. !-----------------------------------------------------------------------
  1748. if (all(weights == zero)) return
  1749. !-----------------------------------------------------------------------
  1750. !
  1751. ! restrict the range of links to search for existing links
  1752. !
  1753. !-----------------------------------------------------------------------
  1754. if (first_call) then
  1755. allocate(link_add1(2,grid1_size), link_add2(2,grid2_size))
  1756. link_add1 = 0
  1757. link_add2 = 0
  1758. first_call = .false.
  1759. min_link = 1
  1760. max_link = 0
  1761. else
  1762. min_link = min(link_add1(1,add1),link_add2(1,add2))
  1763. max_link = max(link_add1(2,add1),link_add2(2,add2))
  1764. if (min_link == 0) then
  1765. min_link = 1
  1766. max_link = 0
  1767. endif
  1768. endif
  1769. !-----------------------------------------------------------------------
  1770. !
  1771. ! if the link already exists, add the weight to the current weight
  1772. ! arrays
  1773. !
  1774. !-----------------------------------------------------------------------
  1775. do nlink=min_link,max_link
  1776. if (add1 == grid1_add_map1(nlink)) then
  1777. if (add2 == grid2_add_map1(nlink)) then
  1778. wts_map1(:,nlink) = wts_map1(:,nlink) + weights(1:num_wts)
  1779. if (num_maps == 2) then
  1780. wts_map2(:,nlink) = wts_map2(:,nlink) +
  1781. & weights(num_wts+1:2*num_wts)
  1782. endif
  1783. return
  1784. endif
  1785. endif
  1786. end do
  1787. !-----------------------------------------------------------------------
  1788. !
  1789. ! if the link does not yet exist, increment number of links and
  1790. ! check to see if remap arrays need to be increased to accomodate
  1791. ! the new link. then store the link.
  1792. !
  1793. !-----------------------------------------------------------------------
  1794. num_links_map1 = num_links_map1 + 1
  1795. if (num_links_map1 > max_links_map1)
  1796. & call resize_remap_vars(1,resize_increment)
  1797. grid1_add_map1(num_links_map1) = add1
  1798. grid2_add_map1(num_links_map1) = add2
  1799. wts_map1 (:,num_links_map1) = weights(1:num_wts)
  1800. if (num_maps > 1) then
  1801. num_links_map2 = num_links_map2 + 1
  1802. if (num_links_map2 > max_links_map2)
  1803. & call resize_remap_vars(2,resize_increment)
  1804. grid1_add_map2(num_links_map2) = add1
  1805. grid2_add_map2(num_links_map2) = add2
  1806. wts_map2 (:,num_links_map2) = weights(num_wts+1:2*num_wts)
  1807. endif
  1808. if (link_add1(1,add1) == 0) link_add1(1,add1) = num_links_map1
  1809. if (link_add2(1,add2) == 0) link_add2(1,add2) = num_links_map1
  1810. link_add1(2,add1) = num_links_map1
  1811. link_add2(2,add2) = num_links_map1
  1812. !-----------------------------------------------------------------------
  1813. end subroutine store_link_cnsrv
  1814. !***********************************************************************
  1815. end module remap_conservative
  1816. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!