remap_write.f 56 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475147614771478147914801481148214831484148514861487148814891490149114921493149414951496149714981499150015011502150315041505150615071508150915101511151215131514151515161517151815191520152115221523152415251526152715281529153015311532153315341535153615371538153915401541154215431544154515461547154815491550155115521553155415551556155715581559156015611562156315641565156615671568156915701571157215731574157515761577157815791580158115821583158415851586158715881589159015911592159315941595159615971598159916001601160216031604160516061607160816091610161116121613161416151616161716181619162016211622162316241625162616271628162916301631163216331634163516361637163816391640164116421643164416451646164716481649165016511652165316541655165616571658165916601661166216631664166516661667166816691670167116721673167416751676167716781679168016811682168316841685168616871688168916901691169216931694169516961697169816991700170117021703170417051706170717081709171017111712171317141715171617171718171917201721172217231724172517261727172817291730173117321733173417351736173717381739174017411742174317441745174617471748174917501751175217531754175517561757175817591760176117621763
  1. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  2. !
  3. ! This module contains routines for writing the remapping data to
  4. ! a file. Before writing the data for each mapping, the links are
  5. ! sorted by destination grid address.
  6. !
  7. !-----------------------------------------------------------------------
  8. !
  9. ! CVS:$Id: remap_write.f,v 1.7 2001/08/21 21:06:42 pwjones Exp $
  10. !
  11. ! Copyright (c) 1997, 1998 the Regents of the University of
  12. ! California.
  13. !
  14. ! This software and ancillary information (herein called software)
  15. ! called SCRIP is made available under the terms described here.
  16. ! The software has been approved for release with associated
  17. ! LA-CC Number 98-45.
  18. !
  19. ! Unless otherwise indicated, this software has been authored
  20. ! by an employee or employees of the University of California,
  21. ! operator of the Los Alamos National Laboratory under Contract
  22. ! No. W-7405-ENG-36 with the U.S. Department of Energy. The U.S.
  23. ! Government has rights to use, reproduce, and distribute this
  24. ! software. The public may copy and use this software without
  25. ! charge, provided that this Notice and any statement of authorship
  26. ! are reproduced on all copies. Neither the Government nor the
  27. ! University makes any warranty, express or implied, or assumes
  28. ! any liability or responsibility for the use of this software.
  29. !
  30. ! If software is modified to produce derivative works, such modified
  31. ! software should be clearly marked, so as not to confuse it with
  32. ! the version available from Los Alamos National Laboratory.
  33. !
  34. !***********************************************************************
  35. module remap_write
  36. !-----------------------------------------------------------------------
  37. use kinds_mod ! defines common data types
  38. use constants ! defines common scalar constants
  39. use grids ! module containing grid information
  40. use remap_vars ! module containing remap information
  41. use netcdf_mod ! module with netCDF stuff
  42. implicit none
  43. !-----------------------------------------------------------------------
  44. !
  45. ! module variables
  46. !
  47. !-----------------------------------------------------------------------
  48. character(char_len), private ::
  49. & map_method ! character string for map_type
  50. &, normalize_opt ! character string for normalization option
  51. &, history ! character string for history information
  52. &, convention ! character string for output convention
  53. character(8), private ::
  54. & cdate ! character date string
  55. integer (kind=int_kind), dimension(:), allocatable, private ::
  56. & src_mask_int ! integer masks to determine
  57. &, dst_mask_int ! cells that participate in map
  58. !-----------------------------------------------------------------------
  59. !
  60. ! various netCDF identifiers used by output routines
  61. !
  62. !-----------------------------------------------------------------------
  63. integer (kind=int_kind), private ::
  64. & ncstat ! error flag for netCDF calls
  65. &, nc_file_id ! id for netCDF file
  66. &, nc_srcgrdsize_id ! id for source grid size
  67. &, nc_dstgrdsize_id ! id for destination grid size
  68. &, nc_srcgrdcorn_id ! id for number of source grid corners
  69. &, nc_dstgrdcorn_id ! id for number of dest grid corners
  70. &, nc_srcgrdrank_id ! id for source grid rank
  71. &, nc_dstgrdrank_id ! id for dest grid rank
  72. &, nc_numlinks_id ! id for number of links in mapping
  73. &, nc_numwgts_id ! id for number of weights for mapping
  74. &, nc_srcgrddims_id ! id for source grid dimensions
  75. &, nc_dstgrddims_id ! id for dest grid dimensions
  76. &, nc_srcgrdcntrlat_id ! id for source grid center latitude
  77. &, nc_dstgrdcntrlat_id ! id for dest grid center latitude
  78. &, nc_srcgrdcntrlon_id ! id for source grid center longitude
  79. &, nc_dstgrdcntrlon_id ! id for dest grid center longitude
  80. &, nc_srcgrdimask_id ! id for source grid mask
  81. &, nc_dstgrdimask_id ! id for dest grid mask
  82. &, nc_srcgrdcrnrlat_id ! id for latitude of source grid corners
  83. &, nc_srcgrdcrnrlon_id ! id for longitude of source grid corners
  84. &, nc_dstgrdcrnrlat_id ! id for latitude of dest grid corners
  85. &, nc_dstgrdcrnrlon_id ! id for longitude of dest grid corners
  86. &, nc_srcgrdarea_id ! id for area of source grid cells
  87. &, nc_dstgrdarea_id ! id for area of dest grid cells
  88. &, nc_srcgrdfrac_id ! id for area fraction on source grid
  89. &, nc_dstgrdfrac_id ! id for area fraction on dest grid
  90. &, nc_srcadd_id ! id for map source address
  91. &, nc_dstadd_id ! id for map destination address
  92. &, nc_rmpmatrix_id ! id for remapping matrix
  93. integer (kind=int_kind), dimension(2), private ::
  94. & nc_dims2_id ! netCDF ids for 2d array dims
  95. !***********************************************************************
  96. contains
  97. !***********************************************************************
  98. subroutine write_remap(map1_name, map2_name,
  99. & interp_file1, interp_file2, output_opt)
  100. !-----------------------------------------------------------------------
  101. !
  102. ! calls correct output routine based on output format choice
  103. !
  104. !-----------------------------------------------------------------------
  105. !-----------------------------------------------------------------------
  106. !
  107. ! input variables
  108. !
  109. !-----------------------------------------------------------------------
  110. character(char_len), intent(in) ::
  111. & map1_name, ! name for mapping grid1 to grid2
  112. & map2_name, ! name for mapping grid2 to grid1
  113. & interp_file1, ! filename for map1 remap data
  114. & interp_file2, ! filename for map2 remap data
  115. & output_opt ! option for output conventions
  116. !-----------------------------------------------------------------------
  117. !
  118. ! local variables
  119. !
  120. !-----------------------------------------------------------------------
  121. !-----------------------------------------------------------------------
  122. !
  123. ! define some common variables to be used in all routines
  124. !
  125. !-----------------------------------------------------------------------
  126. select case(norm_opt)
  127. case (norm_opt_none)
  128. normalize_opt = 'none'
  129. case (norm_opt_frcarea)
  130. normalize_opt = 'fracarea'
  131. case (norm_opt_dstarea)
  132. normalize_opt = 'destarea'
  133. end select
  134. select case(map_type)
  135. case(map_type_conserv)
  136. map_method = 'Conservative remapping'
  137. case(map_type_bilinear)
  138. map_method = 'Bilinear remapping'
  139. case(map_type_distwgt)
  140. map_method = 'Distance weighted avg of nearest neighbors'
  141. case(map_type_bicubic)
  142. map_method = 'Bicubic remapping'
  143. case default
  144. stop 'Invalid Map Type'
  145. end select
  146. call date_and_time(date=cdate)
  147. write (history,1000) cdate(5:6),cdate(7:8),cdate(1:4)
  148. 1000 format('Created: ',a2,'-',a2,'-',a4)
  149. !-----------------------------------------------------------------------
  150. !
  151. ! sort address and weight arrays
  152. !
  153. !-----------------------------------------------------------------------
  154. call sort_add(grid2_add_map1, grid1_add_map1, wts_map1)
  155. if (num_maps > 1) then
  156. call sort_add(grid1_add_map2, grid2_add_map2, wts_map2)
  157. endif
  158. !-----------------------------------------------------------------------
  159. !
  160. ! call appropriate output routine
  161. !
  162. !-----------------------------------------------------------------------
  163. select case(output_opt)
  164. case ('scrip')
  165. call write_remap_scrip(map1_name, interp_file1, 1)
  166. case ('ncar-csm')
  167. call write_remap_csm (map1_name, interp_file1, 1)
  168. case default
  169. stop 'unknown output file convention'
  170. end select
  171. !-----------------------------------------------------------------------
  172. !
  173. ! call appropriate output routine for second mapping if required
  174. !
  175. !-----------------------------------------------------------------------
  176. if (num_maps > 1) then
  177. select case(output_opt)
  178. case ('scrip')
  179. call write_remap_scrip(map2_name, interp_file2, 2)
  180. case ('ncar-csm')
  181. call write_remap_csm (map2_name, interp_file2, 2)
  182. case default
  183. stop 'unknown output file convention'
  184. end select
  185. endif
  186. !-----------------------------------------------------------------------
  187. end subroutine write_remap
  188. !***********************************************************************
  189. subroutine write_remap_scrip(map_name, interp_file, direction)
  190. !-----------------------------------------------------------------------
  191. !
  192. ! writes remap data to a netCDF file using SCRIP conventions
  193. !
  194. !-----------------------------------------------------------------------
  195. !-----------------------------------------------------------------------
  196. !
  197. ! input variables
  198. !
  199. !-----------------------------------------------------------------------
  200. character(char_len), intent(in) ::
  201. & map_name ! name for mapping
  202. &, interp_file ! filename for remap data
  203. integer (kind=int_kind), intent(in) ::
  204. & direction ! direction of map (1=grid1 to grid2
  205. ! 2=grid2 to grid1)
  206. !-----------------------------------------------------------------------
  207. !
  208. ! local variables
  209. !
  210. !-----------------------------------------------------------------------
  211. character(char_len) ::
  212. & grid1_ctmp ! character temp for grid1 names
  213. &, grid2_ctmp ! character temp for grid2 names
  214. integer (kind=int_kind) ::
  215. & itmp1 ! integer temp
  216. &, itmp2 ! integer temp
  217. &, itmp3 ! integer temp
  218. &, itmp4 ! integer temp
  219. !-----------------------------------------------------------------------
  220. !
  221. ! create netCDF file for mapping and define some global attributes
  222. !
  223. !-----------------------------------------------------------------------
  224. ncstat = nf_create (interp_file, NF_CLOBBER, nc_file_id)
  225. call netcdf_error_handler(ncstat)
  226. !***
  227. !*** map name
  228. !***
  229. ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, 'title',
  230. & len_trim(map_name), map_name)
  231. call netcdf_error_handler(ncstat)
  232. !***
  233. !*** normalization option
  234. !***
  235. ncstat = nf_put_att_text(nc_file_id, NF_GLOBAL, 'normalization',
  236. & len_trim(normalize_opt), normalize_opt)
  237. call netcdf_error_handler(ncstat)
  238. !***
  239. !*** map method
  240. !***
  241. ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, 'map_method',
  242. & len_trim(map_method), map_method)
  243. call netcdf_error_handler(ncstat)
  244. !***
  245. !*** history
  246. !***
  247. ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, 'history',
  248. & len_trim(history), history)
  249. call netcdf_error_handler(ncstat)
  250. !***
  251. !*** file convention
  252. !***
  253. convention = 'SCRIP'
  254. ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, 'conventions',
  255. & len_trim(convention), convention)
  256. call netcdf_error_handler(ncstat)
  257. !***
  258. !*** source and destination grid names
  259. !***
  260. if (direction == 1) then
  261. grid1_ctmp = 'source_grid'
  262. grid2_ctmp = 'dest_grid'
  263. else
  264. grid1_ctmp = 'dest_grid'
  265. grid2_ctmp = 'source_grid'
  266. endif
  267. ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, trim(grid1_ctmp),
  268. & len_trim(grid1_name), grid1_name)
  269. call netcdf_error_handler(ncstat)
  270. ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, trim(grid2_ctmp),
  271. & len_trim(grid2_name), grid2_name)
  272. call netcdf_error_handler(ncstat)
  273. !-----------------------------------------------------------------------
  274. !
  275. ! prepare netCDF dimension info
  276. !
  277. !-----------------------------------------------------------------------
  278. !***
  279. !*** define grid size dimensions
  280. !***
  281. if (direction == 1) then
  282. itmp1 = grid1_size
  283. itmp2 = grid2_size
  284. else
  285. itmp1 = grid2_size
  286. itmp2 = grid1_size
  287. endif
  288. ncstat = nf_def_dim (nc_file_id, 'src_grid_size', itmp1,
  289. & nc_srcgrdsize_id)
  290. call netcdf_error_handler(ncstat)
  291. ncstat = nf_def_dim (nc_file_id, 'dst_grid_size', itmp2,
  292. & nc_dstgrdsize_id)
  293. call netcdf_error_handler(ncstat)
  294. !***
  295. !*** define grid corner dimension
  296. !***
  297. if (direction == 1) then
  298. itmp1 = grid1_corners
  299. itmp2 = grid2_corners
  300. else
  301. itmp1 = grid2_corners
  302. itmp2 = grid1_corners
  303. endif
  304. ncstat = nf_def_dim (nc_file_id, 'src_grid_corners',
  305. & itmp1, nc_srcgrdcorn_id)
  306. call netcdf_error_handler(ncstat)
  307. ncstat = nf_def_dim (nc_file_id, 'dst_grid_corners',
  308. & itmp2, nc_dstgrdcorn_id)
  309. call netcdf_error_handler(ncstat)
  310. !***
  311. !*** define grid rank dimension
  312. !***
  313. if (direction == 1) then
  314. itmp1 = grid1_rank
  315. itmp2 = grid2_rank
  316. else
  317. itmp1 = grid2_rank
  318. itmp2 = grid1_rank
  319. endif
  320. ncstat = nf_def_dim (nc_file_id, 'src_grid_rank',
  321. & itmp1, nc_srcgrdrank_id)
  322. call netcdf_error_handler(ncstat)
  323. ncstat = nf_def_dim (nc_file_id, 'dst_grid_rank',
  324. & itmp2, nc_dstgrdrank_id)
  325. call netcdf_error_handler(ncstat)
  326. !***
  327. !*** define map size dimensions
  328. !***
  329. if (direction == 1) then
  330. itmp1 = num_links_map1
  331. else
  332. itmp1 = num_links_map2
  333. endif
  334. ncstat = nf_def_dim (nc_file_id, 'num_links',
  335. & itmp1, nc_numlinks_id)
  336. call netcdf_error_handler(ncstat)
  337. ncstat = nf_def_dim (nc_file_id, 'num_wgts',
  338. & num_wts, nc_numwgts_id)
  339. call netcdf_error_handler(ncstat)
  340. !***
  341. !*** define grid dimensions
  342. !***
  343. ncstat = nf_def_var (nc_file_id, 'src_grid_dims', NF_INT,
  344. & 1, nc_srcgrdrank_id, nc_srcgrddims_id)
  345. call netcdf_error_handler(ncstat)
  346. ncstat = nf_def_var (nc_file_id, 'dst_grid_dims', NF_INT,
  347. & 1, nc_dstgrdrank_id, nc_dstgrddims_id)
  348. call netcdf_error_handler(ncstat)
  349. !-----------------------------------------------------------------------
  350. !
  351. ! define all arrays for netCDF descriptors
  352. !
  353. !-----------------------------------------------------------------------
  354. !***
  355. !*** define grid center latitude array
  356. !***
  357. ncstat = nf_def_var (nc_file_id, 'src_grid_center_lat',
  358. & NF_DOUBLE, 1, nc_srcgrdsize_id,
  359. & nc_srcgrdcntrlat_id)
  360. call netcdf_error_handler(ncstat)
  361. ncstat = nf_def_var (nc_file_id, 'dst_grid_center_lat',
  362. & NF_DOUBLE, 1, nc_dstgrdsize_id,
  363. & nc_dstgrdcntrlat_id)
  364. call netcdf_error_handler(ncstat)
  365. !***
  366. !*** define grid center longitude array
  367. !***
  368. ncstat = nf_def_var (nc_file_id, 'src_grid_center_lon',
  369. & NF_DOUBLE, 1, nc_srcgrdsize_id,
  370. & nc_srcgrdcntrlon_id)
  371. call netcdf_error_handler(ncstat)
  372. ncstat = nf_def_var (nc_file_id, 'dst_grid_center_lon',
  373. & NF_DOUBLE, 1, nc_dstgrdsize_id,
  374. & nc_dstgrdcntrlon_id)
  375. call netcdf_error_handler(ncstat)
  376. !***
  377. !*** define grid corner lat/lon arrays
  378. !***
  379. nc_dims2_id(1) = nc_srcgrdcorn_id
  380. nc_dims2_id(2) = nc_srcgrdsize_id
  381. ncstat = nf_def_var (nc_file_id, 'src_grid_corner_lat',
  382. & NF_DOUBLE, 2, nc_dims2_id,
  383. & nc_srcgrdcrnrlat_id)
  384. call netcdf_error_handler(ncstat)
  385. ncstat = nf_def_var (nc_file_id, 'src_grid_corner_lon',
  386. & NF_DOUBLE, 2, nc_dims2_id,
  387. & nc_srcgrdcrnrlon_id)
  388. call netcdf_error_handler(ncstat)
  389. nc_dims2_id(1) = nc_dstgrdcorn_id
  390. nc_dims2_id(2) = nc_dstgrdsize_id
  391. ncstat = nf_def_var (nc_file_id, 'dst_grid_corner_lat',
  392. & NF_DOUBLE, 2, nc_dims2_id,
  393. & nc_dstgrdcrnrlat_id)
  394. call netcdf_error_handler(ncstat)
  395. ncstat = nf_def_var (nc_file_id, 'dst_grid_corner_lon',
  396. & NF_DOUBLE, 2, nc_dims2_id,
  397. & nc_dstgrdcrnrlon_id)
  398. call netcdf_error_handler(ncstat)
  399. !***
  400. !*** define units for all coordinate arrays
  401. !***
  402. if (direction == 1) then
  403. grid1_ctmp = grid1_units
  404. grid2_ctmp = grid2_units
  405. else
  406. grid1_ctmp = grid2_units
  407. grid2_ctmp = grid1_units
  408. endif
  409. ncstat = nf_put_att_text (nc_file_id, nc_srcgrdcntrlat_id,
  410. & 'units', 7, grid1_ctmp)
  411. call netcdf_error_handler(ncstat)
  412. ncstat = nf_put_att_text (nc_file_id, nc_dstgrdcntrlat_id,
  413. & 'units', 7, grid2_ctmp)
  414. call netcdf_error_handler(ncstat)
  415. ncstat = nf_put_att_text (nc_file_id, nc_srcgrdcntrlon_id,
  416. & 'units', 7, grid1_ctmp)
  417. call netcdf_error_handler(ncstat)
  418. ncstat = nf_put_att_text (nc_file_id, nc_dstgrdcntrlon_id,
  419. & 'units', 7, grid2_ctmp)
  420. call netcdf_error_handler(ncstat)
  421. ncstat = nf_put_att_text (nc_file_id, nc_srcgrdcrnrlat_id,
  422. & 'units', 7, grid1_ctmp)
  423. call netcdf_error_handler(ncstat)
  424. ncstat = nf_put_att_text (nc_file_id, nc_srcgrdcrnrlon_id,
  425. & 'units', 7, grid1_ctmp)
  426. call netcdf_error_handler(ncstat)
  427. ncstat = nf_put_att_text (nc_file_id, nc_dstgrdcrnrlat_id,
  428. & 'units', 7, grid2_ctmp)
  429. call netcdf_error_handler(ncstat)
  430. ncstat = nf_put_att_text (nc_file_id, nc_dstgrdcrnrlon_id,
  431. & 'units', 7, grid2_ctmp)
  432. call netcdf_error_handler(ncstat)
  433. !***
  434. !*** define grid mask
  435. !***
  436. ncstat = nf_def_var (nc_file_id, 'src_grid_imask', NF_INT,
  437. & 1, nc_srcgrdsize_id, nc_srcgrdimask_id)
  438. call netcdf_error_handler(ncstat)
  439. ncstat = nf_put_att_text (nc_file_id, nc_srcgrdimask_id,
  440. & 'units', 8, 'unitless')
  441. call netcdf_error_handler(ncstat)
  442. ncstat = nf_def_var (nc_file_id, 'dst_grid_imask', NF_INT,
  443. & 1, nc_dstgrdsize_id, nc_dstgrdimask_id)
  444. call netcdf_error_handler(ncstat)
  445. ncstat = nf_put_att_text (nc_file_id, nc_dstgrdimask_id,
  446. & 'units', 8, 'unitless')
  447. call netcdf_error_handler(ncstat)
  448. !***
  449. !*** define grid area arrays
  450. !***
  451. ncstat = nf_def_var (nc_file_id, 'src_grid_area',
  452. & NF_DOUBLE, 1, nc_srcgrdsize_id,
  453. & nc_srcgrdarea_id)
  454. call netcdf_error_handler(ncstat)
  455. ncstat = nf_put_att_text (nc_file_id, nc_srcgrdarea_id,
  456. & 'units', 14, 'square radians')
  457. call netcdf_error_handler(ncstat)
  458. ncstat = nf_def_var (nc_file_id, 'dst_grid_area',
  459. & NF_DOUBLE, 1, nc_dstgrdsize_id,
  460. & nc_dstgrdarea_id)
  461. call netcdf_error_handler(ncstat)
  462. ncstat = nf_put_att_text (nc_file_id, nc_dstgrdarea_id,
  463. & 'units', 14, 'square radians')
  464. call netcdf_error_handler(ncstat)
  465. !***
  466. !*** define grid fraction arrays
  467. !***
  468. ncstat = nf_def_var (nc_file_id, 'src_grid_frac',
  469. & NF_DOUBLE, 1, nc_srcgrdsize_id,
  470. & nc_srcgrdfrac_id)
  471. call netcdf_error_handler(ncstat)
  472. ncstat = nf_put_att_text (nc_file_id, nc_srcgrdfrac_id,
  473. & 'units', 8, 'unitless')
  474. call netcdf_error_handler(ncstat)
  475. ncstat = nf_def_var (nc_file_id, 'dst_grid_frac',
  476. & NF_DOUBLE, 1, nc_dstgrdsize_id,
  477. & nc_dstgrdfrac_id)
  478. call netcdf_error_handler(ncstat)
  479. ncstat = nf_put_att_text (nc_file_id, nc_dstgrdfrac_id,
  480. & 'units', 8, 'unitless')
  481. call netcdf_error_handler(ncstat)
  482. !***
  483. !*** define mapping arrays
  484. !***
  485. ncstat = nf_def_var (nc_file_id, 'src_address',
  486. & NF_INT, 1, nc_numlinks_id,
  487. & nc_srcadd_id)
  488. call netcdf_error_handler(ncstat)
  489. ncstat = nf_def_var (nc_file_id, 'dst_address',
  490. & NF_INT, 1, nc_numlinks_id,
  491. & nc_dstadd_id)
  492. call netcdf_error_handler(ncstat)
  493. nc_dims2_id(1) = nc_numwgts_id
  494. nc_dims2_id(2) = nc_numlinks_id
  495. ncstat = nf_def_var (nc_file_id, 'remap_matrix',
  496. & NF_DOUBLE, 2, nc_dims2_id,
  497. & nc_rmpmatrix_id)
  498. call netcdf_error_handler(ncstat)
  499. !***
  500. !*** end definition stage
  501. !***
  502. ncstat = nf_enddef(nc_file_id)
  503. call netcdf_error_handler(ncstat)
  504. !-----------------------------------------------------------------------
  505. !
  506. ! compute integer masks
  507. !
  508. !-----------------------------------------------------------------------
  509. if (direction == 1) then
  510. allocate (src_mask_int(grid1_size),
  511. & dst_mask_int(grid2_size))
  512. where (grid2_mask)
  513. dst_mask_int = 1
  514. elsewhere
  515. dst_mask_int = 0
  516. endwhere
  517. where (grid1_mask)
  518. src_mask_int = 1
  519. elsewhere
  520. src_mask_int = 0
  521. endwhere
  522. else
  523. allocate (src_mask_int(grid2_size),
  524. & dst_mask_int(grid1_size))
  525. where (grid1_mask)
  526. dst_mask_int = 1
  527. elsewhere
  528. dst_mask_int = 0
  529. endwhere
  530. where (grid2_mask)
  531. src_mask_int = 1
  532. elsewhere
  533. src_mask_int = 0
  534. endwhere
  535. endif
  536. !-----------------------------------------------------------------------
  537. !
  538. ! change units of lat/lon coordinates if input units different
  539. ! from radians
  540. !
  541. !-----------------------------------------------------------------------
  542. if (grid1_units(1:7) == 'degrees' .and. direction == 1) then
  543. grid1_center_lat = grid1_center_lat/deg2rad
  544. grid1_center_lon = grid1_center_lon/deg2rad
  545. grid1_corner_lat = grid1_corner_lat/deg2rad
  546. grid1_corner_lon = grid1_corner_lon/deg2rad
  547. endif
  548. if (grid2_units(1:7) == 'degrees' .and. direction == 1) then
  549. grid2_center_lat = grid2_center_lat/deg2rad
  550. grid2_center_lon = grid2_center_lon/deg2rad
  551. grid2_corner_lat = grid2_corner_lat/deg2rad
  552. grid2_corner_lon = grid2_corner_lon/deg2rad
  553. endif
  554. !-----------------------------------------------------------------------
  555. !
  556. ! write mapping data
  557. !
  558. !-----------------------------------------------------------------------
  559. if (direction == 1) then
  560. itmp1 = nc_srcgrddims_id
  561. itmp2 = nc_dstgrddims_id
  562. else
  563. itmp2 = nc_srcgrddims_id
  564. itmp1 = nc_dstgrddims_id
  565. endif
  566. ncstat = nf_put_var_int(nc_file_id, itmp1, grid1_dims)
  567. call netcdf_error_handler(ncstat)
  568. ncstat = nf_put_var_int(nc_file_id, itmp2, grid2_dims)
  569. call netcdf_error_handler(ncstat)
  570. ncstat = nf_put_var_int(nc_file_id, nc_srcgrdimask_id,
  571. & src_mask_int)
  572. call netcdf_error_handler(ncstat)
  573. ncstat = nf_put_var_int(nc_file_id, nc_dstgrdimask_id,
  574. & dst_mask_int)
  575. call netcdf_error_handler(ncstat)
  576. deallocate(src_mask_int, dst_mask_int)
  577. if (direction == 1) then
  578. itmp1 = nc_srcgrdcntrlat_id
  579. itmp2 = nc_srcgrdcntrlon_id
  580. itmp3 = nc_srcgrdcrnrlat_id
  581. itmp4 = nc_srcgrdcrnrlon_id
  582. else
  583. itmp1 = nc_dstgrdcntrlat_id
  584. itmp2 = nc_dstgrdcntrlon_id
  585. itmp3 = nc_dstgrdcrnrlat_id
  586. itmp4 = nc_dstgrdcrnrlon_id
  587. endif
  588. ncstat = nf_put_var_double(nc_file_id, itmp1, grid1_center_lat)
  589. call netcdf_error_handler(ncstat)
  590. ncstat = nf_put_var_double(nc_file_id, itmp2, grid1_center_lon)
  591. call netcdf_error_handler(ncstat)
  592. ncstat = nf_put_var_double(nc_file_id, itmp3, grid1_corner_lat)
  593. call netcdf_error_handler(ncstat)
  594. ncstat = nf_put_var_double(nc_file_id, itmp4, grid1_corner_lon)
  595. call netcdf_error_handler(ncstat)
  596. if (direction == 1) then
  597. itmp1 = nc_dstgrdcntrlat_id
  598. itmp2 = nc_dstgrdcntrlon_id
  599. itmp3 = nc_dstgrdcrnrlat_id
  600. itmp4 = nc_dstgrdcrnrlon_id
  601. else
  602. itmp1 = nc_srcgrdcntrlat_id
  603. itmp2 = nc_srcgrdcntrlon_id
  604. itmp3 = nc_srcgrdcrnrlat_id
  605. itmp4 = nc_srcgrdcrnrlon_id
  606. endif
  607. ncstat = nf_put_var_double(nc_file_id, itmp1, grid2_center_lat)
  608. call netcdf_error_handler(ncstat)
  609. ncstat = nf_put_var_double(nc_file_id, itmp2, grid2_center_lon)
  610. call netcdf_error_handler(ncstat)
  611. ncstat = nf_put_var_double(nc_file_id, itmp3, grid2_corner_lat)
  612. call netcdf_error_handler(ncstat)
  613. ncstat = nf_put_var_double(nc_file_id, itmp4, grid2_corner_lon)
  614. call netcdf_error_handler(ncstat)
  615. if (direction == 1) then
  616. itmp1 = nc_srcgrdarea_id
  617. itmp2 = nc_srcgrdfrac_id
  618. itmp3 = nc_dstgrdarea_id
  619. itmp4 = nc_dstgrdfrac_id
  620. else
  621. itmp1 = nc_dstgrdarea_id
  622. itmp2 = nc_dstgrdfrac_id
  623. itmp3 = nc_srcgrdarea_id
  624. itmp4 = nc_srcgrdfrac_id
  625. endif
  626. if (luse_grid1_area) then
  627. ncstat = nf_put_var_double(nc_file_id, itmp1, grid1_area_in)
  628. else
  629. ncstat = nf_put_var_double(nc_file_id, itmp1, grid1_area)
  630. endif
  631. call netcdf_error_handler(ncstat)
  632. ncstat = nf_put_var_double(nc_file_id, itmp2, grid1_frac)
  633. call netcdf_error_handler(ncstat)
  634. if (luse_grid2_area) then
  635. ncstat = nf_put_var_double(nc_file_id, itmp3, grid2_area_in)
  636. else
  637. ncstat = nf_put_var_double(nc_file_id, itmp3, grid2_area)
  638. endif
  639. call netcdf_error_handler(ncstat)
  640. ncstat = nf_put_var_double(nc_file_id, itmp4, grid2_frac)
  641. call netcdf_error_handler(ncstat)
  642. if (direction == 1) then
  643. ncstat = nf_put_var_int(nc_file_id, nc_srcadd_id,
  644. & grid1_add_map1)
  645. call netcdf_error_handler(ncstat)
  646. ncstat = nf_put_var_int(nc_file_id, nc_dstadd_id,
  647. & grid2_add_map1)
  648. call netcdf_error_handler(ncstat)
  649. ncstat = nf_put_var_double(nc_file_id, nc_rmpmatrix_id,
  650. & wts_map1)
  651. call netcdf_error_handler(ncstat)
  652. else
  653. ncstat = nf_put_var_int(nc_file_id, nc_srcadd_id,
  654. & grid2_add_map2)
  655. call netcdf_error_handler(ncstat)
  656. ncstat = nf_put_var_int(nc_file_id, nc_dstadd_id,
  657. & grid1_add_map2)
  658. call netcdf_error_handler(ncstat)
  659. ncstat = nf_put_var_double(nc_file_id, nc_rmpmatrix_id,
  660. & wts_map2)
  661. call netcdf_error_handler(ncstat)
  662. endif
  663. ncstat = nf_close(nc_file_id)
  664. call netcdf_error_handler(ncstat)
  665. !-----------------------------------------------------------------------
  666. end subroutine write_remap_scrip
  667. !***********************************************************************
  668. subroutine write_remap_csm(map_name, interp_file, direction)
  669. !-----------------------------------------------------------------------
  670. !
  671. ! writes remap data to a netCDF file using NCAR-CSM conventions
  672. !
  673. !-----------------------------------------------------------------------
  674. !-----------------------------------------------------------------------
  675. !
  676. ! input variables
  677. !
  678. !-----------------------------------------------------------------------
  679. character(char_len), intent(in) ::
  680. & map_name ! name for mapping
  681. &, interp_file ! filename for remap data
  682. integer (kind=int_kind), intent(in) ::
  683. & direction ! direction of map (1=grid1 to grid2
  684. ! 2=grid2 to grid1)
  685. !-----------------------------------------------------------------------
  686. !
  687. ! local variables
  688. !
  689. !-----------------------------------------------------------------------
  690. character(char_len) ::
  691. & grid1_ctmp ! character temp for grid1 names
  692. &, grid2_ctmp ! character temp for grid2 names
  693. integer (kind=int_kind) ::
  694. & itmp1 ! integer temp
  695. &, itmp2 ! integer temp
  696. &, itmp3 ! integer temp
  697. &, itmp4 ! integer temp
  698. &, nc_numwgts1_id ! extra netCDF id for additional weights
  699. &, nc_src_isize_id ! extra netCDF id for ni_a
  700. &, nc_src_jsize_id ! extra netCDF id for nj_a
  701. &, nc_dst_isize_id ! extra netCDF id for ni_b
  702. &, nc_dst_jsize_id ! extra netCDF id for nj_b
  703. &, nc_rmpmatrix2_id ! extra netCDF id for high-order remap matrix
  704. real (kind=dbl_kind), dimension(:),allocatable ::
  705. & wts1 ! CSM wants single array for 1st-order wts
  706. real (kind=dbl_kind), dimension(:,:),allocatable ::
  707. & wts2 ! write remaining weights in different array
  708. !-----------------------------------------------------------------------
  709. !
  710. ! create netCDF file for mapping and define some global attributes
  711. !
  712. !-----------------------------------------------------------------------
  713. ncstat = nf_create (interp_file, NF_CLOBBER, nc_file_id)
  714. call netcdf_error_handler(ncstat)
  715. !***
  716. !*** map name
  717. !***
  718. ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, 'title',
  719. & len_trim(map_name), map_name)
  720. call netcdf_error_handler(ncstat)
  721. !***
  722. !*** normalization option
  723. !***
  724. ncstat = nf_put_att_text(nc_file_id, NF_GLOBAL, 'normalization',
  725. & len_trim(normalize_opt), normalize_opt)
  726. call netcdf_error_handler(ncstat)
  727. !***
  728. !*** map method
  729. !***
  730. ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, 'map_method',
  731. & len_trim(map_method), map_method)
  732. call netcdf_error_handler(ncstat)
  733. !***
  734. !*** history
  735. !***
  736. ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, 'history',
  737. & len_trim(history), history)
  738. call netcdf_error_handler(ncstat)
  739. !***
  740. !*** file convention
  741. !***
  742. convention = 'NCAR-CSM'
  743. ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, 'conventions',
  744. & len_trim(convention), convention)
  745. call netcdf_error_handler(ncstat)
  746. !***
  747. !*** source and destination grid names
  748. !***
  749. if (direction == 1) then
  750. grid1_ctmp = 'domain_a'
  751. grid2_ctmp = 'domain_b'
  752. else
  753. grid1_ctmp = 'domain_b'
  754. grid2_ctmp = 'domain_a'
  755. endif
  756. ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, trim(grid1_ctmp),
  757. & len_trim(grid1_name), grid1_name)
  758. call netcdf_error_handler(ncstat)
  759. ncstat = nf_put_att_text (nc_file_id, NF_GLOBAL, trim(grid2_ctmp),
  760. & len_trim(grid2_name), grid2_name)
  761. call netcdf_error_handler(ncstat)
  762. !-----------------------------------------------------------------------
  763. !
  764. ! prepare netCDF dimension info
  765. !
  766. !-----------------------------------------------------------------------
  767. !***
  768. !*** define grid size dimensions
  769. !***
  770. if (direction == 1) then
  771. itmp1 = grid1_size
  772. itmp2 = grid2_size
  773. else
  774. itmp1 = grid2_size
  775. itmp2 = grid1_size
  776. endif
  777. ncstat = nf_def_dim (nc_file_id, 'n_a', itmp1, nc_srcgrdsize_id)
  778. call netcdf_error_handler(ncstat)
  779. ncstat = nf_def_dim (nc_file_id, 'n_b', itmp2, nc_dstgrdsize_id)
  780. call netcdf_error_handler(ncstat)
  781. !***
  782. !*** define grid corner dimension
  783. !***
  784. if (direction == 1) then
  785. itmp1 = grid1_corners
  786. itmp2 = grid2_corners
  787. else
  788. itmp1 = grid2_corners
  789. itmp2 = grid1_corners
  790. endif
  791. ncstat = nf_def_dim (nc_file_id, 'nv_a', itmp1, nc_srcgrdcorn_id)
  792. call netcdf_error_handler(ncstat)
  793. ncstat = nf_def_dim (nc_file_id, 'nv_b', itmp2, nc_dstgrdcorn_id)
  794. call netcdf_error_handler(ncstat)
  795. !***
  796. !*** define grid rank dimension
  797. !***
  798. if (direction == 1) then
  799. itmp1 = grid1_rank
  800. itmp2 = grid2_rank
  801. else
  802. itmp1 = grid2_rank
  803. itmp2 = grid1_rank
  804. endif
  805. ncstat = nf_def_dim (nc_file_id, 'src_grid_rank',
  806. & itmp1, nc_srcgrdrank_id)
  807. call netcdf_error_handler(ncstat)
  808. ncstat = nf_def_dim (nc_file_id, 'dst_grid_rank',
  809. & itmp2, nc_dstgrdrank_id)
  810. call netcdf_error_handler(ncstat)
  811. !***
  812. !*** define first two dims as if 2-d cartesian domain
  813. !***
  814. if (direction == 1) then
  815. itmp1 = grid1_dims(1)
  816. if (grid1_rank > 1) then
  817. itmp2 = grid1_dims(2)
  818. else
  819. itmp2 = 0
  820. endif
  821. itmp3 = grid2_dims(1)
  822. if (grid2_rank > 1) then
  823. itmp4 = grid2_dims(2)
  824. else
  825. itmp4 = 0
  826. endif
  827. else
  828. itmp1 = grid2_dims(1)
  829. if (grid2_rank > 1) then
  830. itmp2 = grid2_dims(2)
  831. else
  832. itmp2 = 0
  833. endif
  834. itmp3 = grid1_dims(1)
  835. if (grid1_rank > 1) then
  836. itmp4 = grid1_dims(2)
  837. else
  838. itmp4 = 0
  839. endif
  840. endif
  841. ncstat = nf_def_dim (nc_file_id, 'ni_a', itmp1, nc_src_isize_id)
  842. call netcdf_error_handler(ncstat)
  843. ncstat = nf_def_dim (nc_file_id, 'nj_a', itmp2, nc_src_jsize_id)
  844. call netcdf_error_handler(ncstat)
  845. ncstat = nf_def_dim (nc_file_id, 'ni_b', itmp3, nc_dst_isize_id)
  846. call netcdf_error_handler(ncstat)
  847. ncstat = nf_def_dim (nc_file_id, 'nj_b', itmp4, nc_dst_jsize_id)
  848. call netcdf_error_handler(ncstat)
  849. !***
  850. !*** define map size dimensions
  851. !***
  852. if (direction == 1) then
  853. itmp1 = num_links_map1
  854. else
  855. itmp1 = num_links_map2
  856. endif
  857. ncstat = nf_def_dim (nc_file_id, 'n_s', itmp1, nc_numlinks_id)
  858. call netcdf_error_handler(ncstat)
  859. ncstat = nf_def_dim (nc_file_id, 'num_wgts',
  860. & num_wts, nc_numwgts_id)
  861. call netcdf_error_handler(ncstat)
  862. if (num_wts > 1) then
  863. ncstat = nf_def_dim (nc_file_id, 'num_wgts1',
  864. & num_wts-1, nc_numwgts1_id)
  865. call netcdf_error_handler(ncstat)
  866. endif
  867. !***
  868. !*** define grid dimensions
  869. !***
  870. ncstat = nf_def_var (nc_file_id, 'src_grid_dims', NF_INT,
  871. & 1, nc_srcgrdrank_id, nc_srcgrddims_id)
  872. call netcdf_error_handler(ncstat)
  873. ncstat = nf_def_var (nc_file_id, 'dst_grid_dims', NF_INT,
  874. & 1, nc_dstgrdrank_id, nc_dstgrddims_id)
  875. call netcdf_error_handler(ncstat)
  876. !-----------------------------------------------------------------------
  877. !
  878. ! define all arrays for netCDF descriptors
  879. !
  880. !-----------------------------------------------------------------------
  881. !***
  882. !*** define grid center latitude array
  883. !***
  884. ncstat = nf_def_var (nc_file_id, 'yc_a',
  885. & NF_DOUBLE, 1, nc_srcgrdsize_id,
  886. & nc_srcgrdcntrlat_id)
  887. call netcdf_error_handler(ncstat)
  888. ncstat = nf_def_var (nc_file_id, 'yc_b',
  889. & NF_DOUBLE, 1, nc_dstgrdsize_id,
  890. & nc_dstgrdcntrlat_id)
  891. call netcdf_error_handler(ncstat)
  892. !***
  893. !*** define grid center longitude array
  894. !***
  895. ncstat = nf_def_var (nc_file_id, 'xc_a',
  896. & NF_DOUBLE, 1, nc_srcgrdsize_id,
  897. & nc_srcgrdcntrlon_id)
  898. call netcdf_error_handler(ncstat)
  899. ncstat = nf_def_var (nc_file_id, 'xc_b',
  900. & NF_DOUBLE, 1, nc_dstgrdsize_id,
  901. & nc_dstgrdcntrlon_id)
  902. call netcdf_error_handler(ncstat)
  903. !***
  904. !*** define grid corner lat/lon arrays
  905. !***
  906. nc_dims2_id(1) = nc_srcgrdcorn_id
  907. nc_dims2_id(2) = nc_srcgrdsize_id
  908. ncstat = nf_def_var (nc_file_id, 'yv_a',
  909. & NF_DOUBLE, 2, nc_dims2_id,
  910. & nc_srcgrdcrnrlat_id)
  911. call netcdf_error_handler(ncstat)
  912. ncstat = nf_def_var (nc_file_id, 'xv_a',
  913. & NF_DOUBLE, 2, nc_dims2_id,
  914. & nc_srcgrdcrnrlon_id)
  915. call netcdf_error_handler(ncstat)
  916. nc_dims2_id(1) = nc_dstgrdcorn_id
  917. nc_dims2_id(2) = nc_dstgrdsize_id
  918. ncstat = nf_def_var (nc_file_id, 'yv_b',
  919. & NF_DOUBLE, 2, nc_dims2_id,
  920. & nc_dstgrdcrnrlat_id)
  921. call netcdf_error_handler(ncstat)
  922. ncstat = nf_def_var (nc_file_id, 'xv_b',
  923. & NF_DOUBLE, 2, nc_dims2_id,
  924. & nc_dstgrdcrnrlon_id)
  925. call netcdf_error_handler(ncstat)
  926. !***
  927. !*** CSM wants all in degrees
  928. !***
  929. grid1_units = 'degrees'
  930. grid2_units = 'degrees'
  931. if (direction == 1) then
  932. grid1_ctmp = grid1_units
  933. grid2_ctmp = grid2_units
  934. else
  935. grid1_ctmp = grid2_units
  936. grid2_ctmp = grid1_units
  937. endif
  938. ncstat = nf_put_att_text (nc_file_id, nc_srcgrdcntrlat_id,
  939. & 'units', 7, grid1_ctmp)
  940. call netcdf_error_handler(ncstat)
  941. ncstat = nf_put_att_text (nc_file_id, nc_dstgrdcntrlat_id,
  942. & 'units', 7, grid2_ctmp)
  943. call netcdf_error_handler(ncstat)
  944. ncstat = nf_put_att_text (nc_file_id, nc_srcgrdcntrlon_id,
  945. & 'units', 7, grid1_ctmp)
  946. call netcdf_error_handler(ncstat)
  947. ncstat = nf_put_att_text (nc_file_id, nc_dstgrdcntrlon_id,
  948. & 'units', 7, grid2_ctmp)
  949. call netcdf_error_handler(ncstat)
  950. ncstat = nf_put_att_text (nc_file_id, nc_srcgrdcrnrlat_id,
  951. & 'units', 7, grid1_ctmp)
  952. call netcdf_error_handler(ncstat)
  953. ncstat = nf_put_att_text (nc_file_id, nc_srcgrdcrnrlon_id,
  954. & 'units', 7, grid1_ctmp)
  955. call netcdf_error_handler(ncstat)
  956. ncstat = nf_put_att_text (nc_file_id, nc_dstgrdcrnrlat_id,
  957. & 'units', 7, grid2_ctmp)
  958. call netcdf_error_handler(ncstat)
  959. ncstat = nf_put_att_text (nc_file_id, nc_dstgrdcrnrlon_id,
  960. & 'units', 7, grid2_ctmp)
  961. call netcdf_error_handler(ncstat)
  962. !***
  963. !*** define grid mask
  964. !***
  965. ncstat = nf_def_var (nc_file_id, 'mask_a', NF_INT,
  966. & 1, nc_srcgrdsize_id, nc_srcgrdimask_id)
  967. call netcdf_error_handler(ncstat)
  968. ncstat = nf_put_att_text (nc_file_id, nc_srcgrdimask_id,
  969. & 'units', 8, 'unitless')
  970. call netcdf_error_handler(ncstat)
  971. ncstat = nf_def_var (nc_file_id, 'mask_b', NF_INT,
  972. & 1, nc_dstgrdsize_id, nc_dstgrdimask_id)
  973. call netcdf_error_handler(ncstat)
  974. ncstat = nf_put_att_text (nc_file_id, nc_dstgrdimask_id,
  975. & 'units', 8, 'unitless')
  976. call netcdf_error_handler(ncstat)
  977. !***
  978. !*** define grid area arrays
  979. !***
  980. ncstat = nf_def_var (nc_file_id, 'area_a',
  981. & NF_DOUBLE, 1, nc_srcgrdsize_id,
  982. & nc_srcgrdarea_id)
  983. call netcdf_error_handler(ncstat)
  984. ncstat = nf_put_att_text (nc_file_id, nc_srcgrdarea_id,
  985. & 'units', 14, 'square radians')
  986. call netcdf_error_handler(ncstat)
  987. ncstat = nf_def_var (nc_file_id, 'area_b',
  988. & NF_DOUBLE, 1, nc_dstgrdsize_id,
  989. & nc_dstgrdarea_id)
  990. call netcdf_error_handler(ncstat)
  991. ncstat = nf_put_att_text (nc_file_id, nc_dstgrdarea_id,
  992. & 'units', 14, 'square radians')
  993. call netcdf_error_handler(ncstat)
  994. !***
  995. !*** define grid fraction arrays
  996. !***
  997. ncstat = nf_def_var (nc_file_id, 'frac_a',
  998. & NF_DOUBLE, 1, nc_srcgrdsize_id,
  999. & nc_srcgrdfrac_id)
  1000. call netcdf_error_handler(ncstat)
  1001. ncstat = nf_put_att_text (nc_file_id, nc_srcgrdfrac_id,
  1002. & 'units', 8, 'unitless')
  1003. call netcdf_error_handler(ncstat)
  1004. ncstat = nf_def_var (nc_file_id, 'frac_b',
  1005. & NF_DOUBLE, 1, nc_dstgrdsize_id,
  1006. & nc_dstgrdfrac_id)
  1007. call netcdf_error_handler(ncstat)
  1008. ncstat = nf_put_att_text (nc_file_id, nc_dstgrdfrac_id,
  1009. & 'units', 8, 'unitless')
  1010. call netcdf_error_handler(ncstat)
  1011. !***
  1012. !*** define mapping arrays
  1013. !***
  1014. ncstat = nf_def_var (nc_file_id, 'col',
  1015. & NF_INT, 1, nc_numlinks_id,
  1016. & nc_srcadd_id)
  1017. call netcdf_error_handler(ncstat)
  1018. ncstat = nf_def_var (nc_file_id, 'row',
  1019. & NF_INT, 1, nc_numlinks_id,
  1020. & nc_dstadd_id)
  1021. call netcdf_error_handler(ncstat)
  1022. ncstat = nf_def_var (nc_file_id, 'S',
  1023. & NF_DOUBLE, 1, nc_numlinks_id,
  1024. & nc_rmpmatrix_id)
  1025. call netcdf_error_handler(ncstat)
  1026. if (num_wts > 1) then
  1027. nc_dims2_id(1) = nc_numwgts1_id
  1028. nc_dims2_id(2) = nc_numlinks_id
  1029. ncstat = nf_def_var (nc_file_id, 'S2',
  1030. & NF_DOUBLE, 2, nc_dims2_id,
  1031. & nc_rmpmatrix2_id)
  1032. call netcdf_error_handler(ncstat)
  1033. endif
  1034. !***
  1035. !*** end definition stage
  1036. !***
  1037. ncstat = nf_enddef(nc_file_id)
  1038. call netcdf_error_handler(ncstat)
  1039. !-----------------------------------------------------------------------
  1040. !
  1041. ! compute integer masks
  1042. !
  1043. !-----------------------------------------------------------------------
  1044. if (direction == 1) then
  1045. allocate (src_mask_int(grid1_size),
  1046. & dst_mask_int(grid2_size))
  1047. where (grid2_mask)
  1048. dst_mask_int = 1
  1049. elsewhere
  1050. dst_mask_int = 0
  1051. endwhere
  1052. where (grid1_mask)
  1053. src_mask_int = 1
  1054. elsewhere
  1055. src_mask_int = 0
  1056. endwhere
  1057. else
  1058. allocate (src_mask_int(grid2_size),
  1059. & dst_mask_int(grid1_size))
  1060. where (grid1_mask)
  1061. dst_mask_int = 1
  1062. elsewhere
  1063. dst_mask_int = 0
  1064. endwhere
  1065. where (grid2_mask)
  1066. src_mask_int = 1
  1067. elsewhere
  1068. src_mask_int = 0
  1069. endwhere
  1070. endif
  1071. !-----------------------------------------------------------------------
  1072. !
  1073. ! change units of lat/lon coordinates if input units different
  1074. ! from radians. if this is the second mapping, the conversion has
  1075. ! alread been done.
  1076. !
  1077. !-----------------------------------------------------------------------
  1078. if (grid1_units(1:7) == 'degrees' .and. direction == 1) then
  1079. grid1_center_lat = grid1_center_lat/deg2rad
  1080. grid1_center_lon = grid1_center_lon/deg2rad
  1081. grid1_corner_lat = grid1_corner_lat/deg2rad
  1082. grid1_corner_lon = grid1_corner_lon/deg2rad
  1083. endif
  1084. if (grid2_units(1:7) == 'degrees' .and. direction == 1) then
  1085. grid2_center_lat = grid2_center_lat/deg2rad
  1086. grid2_center_lon = grid2_center_lon/deg2rad
  1087. grid2_corner_lat = grid2_corner_lat/deg2rad
  1088. grid2_corner_lon = grid2_corner_lon/deg2rad
  1089. endif
  1090. !-----------------------------------------------------------------------
  1091. !
  1092. ! write mapping data
  1093. !
  1094. !-----------------------------------------------------------------------
  1095. if (direction == 1) then
  1096. itmp1 = nc_srcgrddims_id
  1097. itmp2 = nc_dstgrddims_id
  1098. else
  1099. itmp2 = nc_srcgrddims_id
  1100. itmp1 = nc_dstgrddims_id
  1101. endif
  1102. ncstat = nf_put_var_int(nc_file_id, itmp1, grid1_dims)
  1103. call netcdf_error_handler(ncstat)
  1104. ncstat = nf_put_var_int(nc_file_id, itmp2, grid2_dims)
  1105. call netcdf_error_handler(ncstat)
  1106. ncstat = nf_put_var_int(nc_file_id, nc_srcgrdimask_id,
  1107. & src_mask_int)
  1108. call netcdf_error_handler(ncstat)
  1109. ncstat = nf_put_var_int(nc_file_id, nc_dstgrdimask_id,
  1110. & dst_mask_int)
  1111. call netcdf_error_handler(ncstat)
  1112. deallocate(src_mask_int, dst_mask_int)
  1113. if (direction == 1) then
  1114. itmp1 = nc_srcgrdcntrlat_id
  1115. itmp2 = nc_srcgrdcntrlon_id
  1116. itmp3 = nc_srcgrdcrnrlat_id
  1117. itmp4 = nc_srcgrdcrnrlon_id
  1118. else
  1119. itmp1 = nc_dstgrdcntrlat_id
  1120. itmp2 = nc_dstgrdcntrlon_id
  1121. itmp3 = nc_dstgrdcrnrlat_id
  1122. itmp4 = nc_dstgrdcrnrlon_id
  1123. endif
  1124. ncstat = nf_put_var_double(nc_file_id, itmp1, grid1_center_lat)
  1125. call netcdf_error_handler(ncstat)
  1126. ncstat = nf_put_var_double(nc_file_id, itmp2, grid1_center_lon)
  1127. call netcdf_error_handler(ncstat)
  1128. ncstat = nf_put_var_double(nc_file_id, itmp3, grid1_corner_lat)
  1129. call netcdf_error_handler(ncstat)
  1130. ncstat = nf_put_var_double(nc_file_id, itmp4, grid1_corner_lon)
  1131. call netcdf_error_handler(ncstat)
  1132. if (direction == 1) then
  1133. itmp1 = nc_dstgrdcntrlat_id
  1134. itmp2 = nc_dstgrdcntrlon_id
  1135. itmp3 = nc_dstgrdcrnrlat_id
  1136. itmp4 = nc_dstgrdcrnrlon_id
  1137. else
  1138. itmp1 = nc_srcgrdcntrlat_id
  1139. itmp2 = nc_srcgrdcntrlon_id
  1140. itmp3 = nc_srcgrdcrnrlat_id
  1141. itmp4 = nc_srcgrdcrnrlon_id
  1142. endif
  1143. ncstat = nf_put_var_double(nc_file_id, itmp1, grid2_center_lat)
  1144. call netcdf_error_handler(ncstat)
  1145. ncstat = nf_put_var_double(nc_file_id, itmp2, grid2_center_lon)
  1146. call netcdf_error_handler(ncstat)
  1147. ncstat = nf_put_var_double(nc_file_id, itmp3, grid2_corner_lat)
  1148. call netcdf_error_handler(ncstat)
  1149. ncstat = nf_put_var_double(nc_file_id, itmp4, grid2_corner_lon)
  1150. call netcdf_error_handler(ncstat)
  1151. if (direction == 1) then
  1152. itmp1 = nc_srcgrdarea_id
  1153. itmp2 = nc_srcgrdfrac_id
  1154. itmp3 = nc_dstgrdarea_id
  1155. itmp4 = nc_dstgrdfrac_id
  1156. else
  1157. itmp1 = nc_dstgrdarea_id
  1158. itmp2 = nc_dstgrdfrac_id
  1159. itmp3 = nc_srcgrdarea_id
  1160. itmp4 = nc_srcgrdfrac_id
  1161. endif
  1162. if (luse_grid1_area) then
  1163. ncstat = nf_put_var_double(nc_file_id, itmp1, grid1_area_in)
  1164. else
  1165. ncstat = nf_put_var_double(nc_file_id, itmp1, grid1_area)
  1166. endif
  1167. call netcdf_error_handler(ncstat)
  1168. ncstat = nf_put_var_double(nc_file_id, itmp2, grid1_frac)
  1169. call netcdf_error_handler(ncstat)
  1170. if (luse_grid2_area) then
  1171. ncstat = nf_put_var_double(nc_file_id, itmp3, grid2_area)
  1172. else
  1173. ncstat = nf_put_var_double(nc_file_id, itmp3, grid2_area)
  1174. endif
  1175. call netcdf_error_handler(ncstat)
  1176. ncstat = nf_put_var_double(nc_file_id, itmp4, grid2_frac)
  1177. call netcdf_error_handler(ncstat)
  1178. if (direction == 1) then
  1179. ncstat = nf_put_var_int(nc_file_id, nc_srcadd_id,
  1180. & grid1_add_map1)
  1181. call netcdf_error_handler(ncstat)
  1182. ncstat = nf_put_var_int(nc_file_id, nc_dstadd_id,
  1183. & grid2_add_map1)
  1184. call netcdf_error_handler(ncstat)
  1185. if (num_wts == 1) then
  1186. ncstat = nf_put_var_double(nc_file_id, nc_rmpmatrix_id,
  1187. & wts_map1)
  1188. call netcdf_error_handler(ncstat)
  1189. else
  1190. allocate(wts1(num_links_map1),wts2(num_wts-1,num_links_map1))
  1191. wts1 = wts_map1(1,:)
  1192. wts2 = wts_map1(2:,:)
  1193. ncstat = nf_put_var_double(nc_file_id, nc_rmpmatrix_id,
  1194. & wts1)
  1195. call netcdf_error_handler(ncstat)
  1196. ncstat = nf_put_var_double(nc_file_id, nc_rmpmatrix2_id,
  1197. & wts2)
  1198. call netcdf_error_handler(ncstat)
  1199. deallocate(wts1,wts2)
  1200. endif
  1201. else
  1202. ncstat = nf_put_var_int(nc_file_id, nc_srcadd_id,
  1203. & grid2_add_map2)
  1204. call netcdf_error_handler(ncstat)
  1205. ncstat = nf_put_var_int(nc_file_id, nc_dstadd_id,
  1206. & grid1_add_map2)
  1207. call netcdf_error_handler(ncstat)
  1208. if (num_wts == 1) then
  1209. ncstat = nf_put_var_double(nc_file_id, nc_rmpmatrix_id,
  1210. & wts_map2)
  1211. call netcdf_error_handler(ncstat)
  1212. else
  1213. allocate(wts1(num_links_map2),wts2(num_wts-1,num_links_map2))
  1214. wts1 = wts_map2(1,:)
  1215. wts2 = wts_map2(2:,:)
  1216. ncstat = nf_put_var_double(nc_file_id, nc_rmpmatrix_id,
  1217. & wts1)
  1218. call netcdf_error_handler(ncstat)
  1219. ncstat = nf_put_var_double(nc_file_id, nc_rmpmatrix2_id,
  1220. & wts2)
  1221. call netcdf_error_handler(ncstat)
  1222. deallocate(wts1,wts2)
  1223. endif
  1224. endif
  1225. ncstat = nf_close(nc_file_id)
  1226. call netcdf_error_handler(ncstat)
  1227. !-----------------------------------------------------------------------
  1228. end subroutine write_remap_csm
  1229. !***********************************************************************
  1230. subroutine sort_add(add1, add2, weights)
  1231. !-----------------------------------------------------------------------
  1232. !
  1233. ! this routine sorts address and weight arrays based on the
  1234. ! destination address with the source address as a secondary
  1235. ! sorting criterion. the method is a standard heap sort.
  1236. !
  1237. !-----------------------------------------------------------------------
  1238. use kinds_mod ! defines common data types
  1239. use constants ! defines common scalar constants
  1240. implicit none
  1241. !-----------------------------------------------------------------------
  1242. !
  1243. ! Input and Output arrays
  1244. !
  1245. !-----------------------------------------------------------------------
  1246. integer (kind=int_kind), intent(inout), dimension(:) ::
  1247. & add1, ! destination address array (num_links)
  1248. & add2 ! source address array
  1249. real (kind=dbl_kind), intent(inout), dimension(:,:) ::
  1250. & weights ! remapping weights (num_wts, num_links)
  1251. !-----------------------------------------------------------------------
  1252. !
  1253. ! local variables
  1254. !
  1255. !-----------------------------------------------------------------------
  1256. integer (kind=int_kind) ::
  1257. & num_links, ! num of links for this mapping
  1258. & num_wts, ! num of weights for this mapping
  1259. & add1_tmp, add2_tmp, ! temp for addresses during swap
  1260. & nwgt,
  1261. & lvl, final_lvl, ! level indexes for heap sort levels
  1262. & chk_lvl1, chk_lvl2, max_lvl
  1263. real (kind=dbl_kind), dimension(SIZE(weights,DIM=1)) ::
  1264. & wgttmp ! temp for holding wts during swap
  1265. !-----------------------------------------------------------------------
  1266. !
  1267. ! determine total number of links to sort and number of weights
  1268. !
  1269. !-----------------------------------------------------------------------
  1270. num_links = SIZE(add1)
  1271. num_wts = SIZE(weights, DIM=1)
  1272. !-----------------------------------------------------------------------
  1273. !
  1274. ! start at the lowest level (N/2) of the tree and sift lower
  1275. ! values to the bottom of the tree, promoting the larger numbers
  1276. !
  1277. !-----------------------------------------------------------------------
  1278. do lvl=num_links/2,1,-1
  1279. final_lvl = lvl
  1280. add1_tmp = add1(lvl)
  1281. add2_tmp = add2(lvl)
  1282. wgttmp(:) = weights(:,lvl)
  1283. !***
  1284. !*** loop until proper level is found for this link, or reach
  1285. !*** bottom
  1286. !***
  1287. sift_loop1: do
  1288. !***
  1289. !*** find the largest of the two daughters
  1290. !***
  1291. chk_lvl1 = 2*final_lvl
  1292. chk_lvl2 = 2*final_lvl+1
  1293. if (chk_lvl1 .EQ. num_links) chk_lvl2 = chk_lvl1
  1294. if ((add1(chk_lvl1) > add1(chk_lvl2)) .OR.
  1295. & ((add1(chk_lvl1) == add1(chk_lvl2)) .AND.
  1296. & (add2(chk_lvl1) > add2(chk_lvl2)))) then
  1297. max_lvl = chk_lvl1
  1298. else
  1299. max_lvl = chk_lvl2
  1300. endif
  1301. !***
  1302. !*** if the parent is greater than both daughters,
  1303. !*** the correct level has been found
  1304. !***
  1305. if ((add1_tmp .GT. add1(max_lvl)) .OR.
  1306. & ((add1_tmp .EQ. add1(max_lvl)) .AND.
  1307. & (add2_tmp .GT. add2(max_lvl)))) then
  1308. add1(final_lvl) = add1_tmp
  1309. add2(final_lvl) = add2_tmp
  1310. weights(:,final_lvl) = wgttmp(:)
  1311. exit sift_loop1
  1312. !***
  1313. !*** otherwise, promote the largest daughter and push
  1314. !*** down one level in the tree. if haven't reached
  1315. !*** the end of the tree, repeat the process. otherwise
  1316. !*** store last values and exit the loop
  1317. !***
  1318. else
  1319. add1(final_lvl) = add1(max_lvl)
  1320. add2(final_lvl) = add2(max_lvl)
  1321. weights(:,final_lvl) = weights(:,max_lvl)
  1322. final_lvl = max_lvl
  1323. if (2*final_lvl > num_links) then
  1324. add1(final_lvl) = add1_tmp
  1325. add2(final_lvl) = add2_tmp
  1326. weights(:,final_lvl) = wgttmp(:)
  1327. exit sift_loop1
  1328. endif
  1329. endif
  1330. end do sift_loop1
  1331. end do
  1332. !-----------------------------------------------------------------------
  1333. !
  1334. ! now that the heap has been sorted, strip off the top (largest)
  1335. ! value and promote the values below
  1336. !
  1337. !-----------------------------------------------------------------------
  1338. do lvl=num_links,3,-1
  1339. !***
  1340. !*** move the top value and insert it into the correct place
  1341. !***
  1342. add1_tmp = add1(lvl)
  1343. add1(lvl) = add1(1)
  1344. add2_tmp = add2(lvl)
  1345. add2(lvl) = add2(1)
  1346. wgttmp(:) = weights(:,lvl)
  1347. weights(:,lvl) = weights(:,1)
  1348. !***
  1349. !*** as above this loop sifts the tmp values down until proper
  1350. !*** level is reached
  1351. !***
  1352. final_lvl = 1
  1353. sift_loop2: do
  1354. !***
  1355. !*** find the largest of the two daughters
  1356. !***
  1357. chk_lvl1 = 2*final_lvl
  1358. chk_lvl2 = 2*final_lvl+1
  1359. if (chk_lvl2 >= lvl) chk_lvl2 = chk_lvl1
  1360. if ((add1(chk_lvl1) > add1(chk_lvl2)) .OR.
  1361. & ((add1(chk_lvl1) == add1(chk_lvl2)) .AND.
  1362. & (add2(chk_lvl1) > add2(chk_lvl2)))) then
  1363. max_lvl = chk_lvl1
  1364. else
  1365. max_lvl = chk_lvl2
  1366. endif
  1367. !***
  1368. !*** if the parent is greater than both daughters,
  1369. !*** the correct level has been found
  1370. !***
  1371. if ((add1_tmp > add1(max_lvl)) .OR.
  1372. & ((add1_tmp == add1(max_lvl)) .AND.
  1373. & (add2_tmp > add2(max_lvl)))) then
  1374. add1(final_lvl) = add1_tmp
  1375. add2(final_lvl) = add2_tmp
  1376. weights(:,final_lvl) = wgttmp(:)
  1377. exit sift_loop2
  1378. !***
  1379. !*** otherwise, promote the largest daughter and push
  1380. !*** down one level in the tree. if haven't reached
  1381. !*** the end of the tree, repeat the process. otherwise
  1382. !*** store last values and exit the loop
  1383. !***
  1384. else
  1385. add1(final_lvl) = add1(max_lvl)
  1386. add2(final_lvl) = add2(max_lvl)
  1387. weights(:,final_lvl) = weights(:,max_lvl)
  1388. final_lvl = max_lvl
  1389. if (2*final_lvl >= lvl) then
  1390. add1(final_lvl) = add1_tmp
  1391. add2(final_lvl) = add2_tmp
  1392. weights(:,final_lvl) = wgttmp(:)
  1393. exit sift_loop2
  1394. endif
  1395. endif
  1396. end do sift_loop2
  1397. end do
  1398. !***
  1399. !*** swap the last two entries
  1400. !***
  1401. add1_tmp = add1(2)
  1402. add1(2) = add1(1)
  1403. add1(1) = add1_tmp
  1404. add2_tmp = add2(2)
  1405. add2(2) = add2(1)
  1406. add2(1) = add2_tmp
  1407. wgttmp (:) = weights(:,2)
  1408. weights(:,2) = weights(:,1)
  1409. weights(:,1) = wgttmp (:)
  1410. !-----------------------------------------------------------------------
  1411. end subroutine sort_add
  1412. !***********************************************************************
  1413. end module remap_write
  1414. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!