icbrst.F90 22 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436
  1. MODULE icbrst
  2. !!======================================================================
  3. !! *** MODULE icbrst ***
  4. !! Ocean physics: read and write iceberg restart files
  5. !!======================================================================
  6. !! History : 3.3.1 ! 2010-01 (Martin&Adcroft) Original code
  7. !! - ! 2011-03 (Madec) Part conversion to NEMO form
  8. !! - ! Removal of mapping from another grid
  9. !! - ! 2011-04 (Alderson) Split into separate modules
  10. !! - ! 2011-04 (Alderson) Restore restart routine
  11. !! - ! Currently needs a fixed processor
  12. !! - ! layout between restarts
  13. !!----------------------------------------------------------------------
  14. !!----------------------------------------------------------------------
  15. !! icb_rst_read : read restart file
  16. !! icb_rst_write : write restart file
  17. !!----------------------------------------------------------------------
  18. USE par_oce ! NEMO parameters
  19. USE dom_oce ! NEMO domain
  20. USE in_out_manager ! NEMO IO routines
  21. USE lib_mpp ! NEMO MPI library, lk_mpp in particular
  22. USE netcdf ! netcdf routines for IO
  23. USE icb_oce ! define iceberg arrays
  24. USE icbutl ! iceberg utility routines
  25. IMPLICIT NONE
  26. PRIVATE
  27. PUBLIC icb_rst_read ! routine called in icbini.F90 module
  28. PUBLIC icb_rst_write ! routine called in icbstp.F90 module
  29. INTEGER :: nlonid, nlatid, nxid, nyid, nuvelid, nvvelid
  30. INTEGER :: nmassid, nthicknessid, nwidthid, nlengthid
  31. INTEGER :: nyearid, ndayid
  32. INTEGER :: nscaling_id, nmass_of_bits_id, nheat_density_id, numberid
  33. INTEGER :: nsiceid, nsheatid, ncalvid, ncalvhid, nkountid
  34. INTEGER :: nret, ncid, nc_dim
  35. INTEGER, DIMENSION(3) :: nstrt3, nlngth3
  36. !!----------------------------------------------------------------------
  37. !! NEMO/OPA 3.3 , NEMO Consortium (2011)
  38. !! $Id: icbrst.F90 2422 2015-06-05 12:04:13Z ufla $
  39. !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
  40. !!----------------------------------------------------------------------
  41. CONTAINS
  42. SUBROUTINE icb_rst_read()
  43. !!----------------------------------------------------------------------
  44. !! *** SUBROUTINE icb_rst_read ***
  45. !!
  46. !! ** Purpose : read a iceberg restart file
  47. !! NB: for this version, we just read back in the restart for this processor
  48. !! so we cannot change the processor layout currently with iceberg code
  49. !!----------------------------------------------------------------------
  50. INTEGER :: idim, ivar, iatt
  51. INTEGER :: jn, iunlim_dim, ibergs_in_file
  52. INTEGER :: iclass
  53. INTEGER, DIMENSION(1) :: istrt, ilngth, idata
  54. INTEGER, DIMENSION(2) :: istrt2, ilngth2
  55. INTEGER, DIMENSION(nkounts) :: idata2
  56. REAL(wp), DIMENSION(1) :: zdata ! need 1d array to read in with
  57. ! start and count arrays
  58. LOGICAL :: ll_found_restart
  59. CHARACTER(len=256) :: cl_path
  60. CHARACTER(len=256) :: cl_filename
  61. CHARACTER(len=NF90_MAX_NAME) :: cl_dname
  62. TYPE(iceberg) :: localberg ! NOT a pointer but an actual local variable
  63. TYPE(point) :: localpt ! NOT a pointer but an actual local variable
  64. !!----------------------------------------------------------------------
  65. ! Find a restart file. Assume iceberg restarts in same directory as ocean restarts.
  66. cl_path = TRIM(cn_ocerst_indir)
  67. IF( cl_path(LEN_TRIM(cl_path):) /= '/' ) cl_path = TRIM(cl_path) // '/'
  68. cl_filename = ' '
  69. IF ( lk_mpp ) THEN
  70. cl_filename = ' '
  71. WRITE( cl_filename, '("restart_icebergs_",I4.4,".nc")' ) narea-1
  72. INQUIRE( file=TRIM(cl_path)//TRIM(cl_filename), exist=ll_found_restart )
  73. ELSE
  74. cl_filename = 'restart_icebergs.nc'
  75. INQUIRE( file=TRIM(cl_path)//TRIM(cl_filename), exist=ll_found_restart )
  76. ENDIF
  77. IF ( .NOT. ll_found_restart) THEN ! only do the following if a file was found
  78. CALL ctl_stop('icebergs: no restart file found')
  79. ENDIF
  80. IF (nn_verbose_level >= 0 .AND. lwp) &
  81. WRITE(numout,'(2a)') 'icebergs, read_restart_bergs: found restart file = ',TRIM(cl_path)//TRIM(cl_filename)
  82. nret = NF90_OPEN(TRIM(cl_path)//TRIM(cl_filename), NF90_NOWRITE, ncid)
  83. IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart_bergs: nf_open failed')
  84. nret = nf90_inquire(ncid, idim, ivar, iatt, iunlim_dim)
  85. IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart_bergs: nf_inquire failed')
  86. IF( iunlim_dim .NE. -1) THEN
  87. nret = nf90_inquire_dimension(ncid, iunlim_dim, cl_dname, ibergs_in_file)
  88. IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart_bergs: nf_inq_dimlen failed')
  89. nret = NF90_INQ_VARID(ncid, 'number', numberid)
  90. nret = NF90_INQ_VARID(ncid, 'mass_scaling', nscaling_id)
  91. nret = NF90_INQ_VARID(ncid, 'xi', nxid)
  92. nret = NF90_INQ_VARID(ncid, 'yj', nyid)
  93. nret = NF90_INQ_VARID(ncid, 'lon', nlonid)
  94. nret = NF90_INQ_VARID(ncid, 'lat', nlatid)
  95. nret = NF90_INQ_VARID(ncid, 'uvel', nuvelid)
  96. nret = NF90_INQ_VARID(ncid, 'vvel', nvvelid)
  97. nret = NF90_INQ_VARID(ncid, 'mass', nmassid)
  98. nret = NF90_INQ_VARID(ncid, 'thickness', nthicknessid)
  99. nret = NF90_INQ_VARID(ncid, 'width', nwidthid)
  100. nret = NF90_INQ_VARID(ncid, 'length', nlengthid)
  101. nret = NF90_INQ_VARID(ncid, 'year', nyearid)
  102. nret = NF90_INQ_VARID(ncid, 'day', ndayid)
  103. nret = NF90_INQ_VARID(ncid, 'mass_of_bits', nmass_of_bits_id)
  104. nret = NF90_INQ_VARID(ncid, 'heat_density', nheat_density_id)
  105. ilngth(1) = 1
  106. istrt2(1) = 1
  107. ilngth2(1) = nkounts
  108. ilngth2(2) = 1
  109. DO jn=1, ibergs_in_file
  110. istrt(1) = jn
  111. istrt2(2) = jn
  112. nret = NF90_GET_VAR(ncid, numberid, idata2, istrt2, ilngth2 )
  113. localberg%number(:) = idata2(:)
  114. nret = NF90_GET_VAR(ncid, nscaling_id, zdata, istrt, ilngth )
  115. localberg%mass_scaling = zdata(1)
  116. nret = NF90_GET_VAR(ncid, nlonid, zdata, istrt, ilngth)
  117. localpt%lon = zdata(1)
  118. nret = NF90_GET_VAR(ncid, nlatid, zdata, istrt, ilngth)
  119. localpt%lat = zdata(1)
  120. IF (nn_verbose_level >= 2 .AND. lwp) THEN
  121. WRITE(numout,'(a,i5,a,2f10.4,a,i5)') 'icebergs, read_restart_bergs: berg ',jn,' is at ', &
  122. localpt%lon,localpt%lat,' on PE ',narea-1
  123. ENDIF
  124. nret = NF90_GET_VAR(ncid, nxid, zdata, istrt, ilngth)
  125. localpt%xi = zdata(1)
  126. nret = NF90_GET_VAR(ncid, nyid, zdata, istrt, ilngth)
  127. localpt%yj = zdata(1)
  128. nret = NF90_GET_VAR(ncid, nuvelid, zdata, istrt, ilngth )
  129. localpt%uvel = zdata(1)
  130. nret = NF90_GET_VAR(ncid, nvvelid, zdata, istrt, ilngth )
  131. localpt%vvel = zdata(1)
  132. nret = NF90_GET_VAR(ncid, nmassid, zdata, istrt, ilngth )
  133. localpt%mass = zdata(1)
  134. nret = NF90_GET_VAR(ncid, nthicknessid, zdata, istrt, ilngth )
  135. localpt%thickness = zdata(1)
  136. nret = NF90_GET_VAR(ncid, nwidthid, zdata, istrt, ilngth )
  137. localpt%width = zdata(1)
  138. nret = NF90_GET_VAR(ncid, nlengthid, zdata, istrt, ilngth )
  139. localpt%length = zdata(1)
  140. nret = NF90_GET_VAR(ncid, nyearid, idata, istrt, ilngth )
  141. localpt%year = idata(1)
  142. nret = NF90_GET_VAR(ncid, ndayid, zdata, istrt, ilngth )
  143. localpt%day = zdata(1)
  144. nret = NF90_GET_VAR(ncid, nmass_of_bits_id, zdata, istrt, ilngth )
  145. localpt%mass_of_bits = zdata(1)
  146. nret = NF90_GET_VAR(ncid, nheat_density_id, zdata, istrt, ilngth )
  147. localpt%heat_density = zdata(1)
  148. !
  149. CALL icb_utl_add( localberg, localpt )
  150. END DO
  151. !
  152. ENDIF
  153. nret = NF90_INQ_DIMID( ncid, 'c', nc_dim )
  154. IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart: nf_inq_dimid c failed')
  155. nret = NF90_INQUIRE_DIMENSION( ncid, nc_dim, cl_dname, iclass )
  156. IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart: nf_inquire_dimension failed')
  157. nret = NF90_INQ_VARID(ncid, 'kount' , nkountid)
  158. nret = NF90_INQ_VARID(ncid, 'calving' , ncalvid)
  159. nret = NF90_INQ_VARID(ncid, 'calving_hflx', ncalvhid)
  160. nret = NF90_INQ_VARID(ncid, 'stored_ice' , nsiceid)
  161. nret = NF90_INQ_VARID(ncid, 'stored_heat' , nsheatid)
  162. nstrt3(1) = 1
  163. nstrt3(2) = 1
  164. nlngth3(1) = jpi
  165. nlngth3(2) = jpj
  166. nlngth3(3) = 1
  167. DO jn = 1, iclass
  168. nstrt3(3) = jn
  169. nret = NF90_GET_VAR( ncid, nsiceid , griddata, nstrt3, nlngth3 )
  170. berg_grid%stored_ice(:,:,jn) = griddata(:,:,1)
  171. END DO
  172. nret = NF90_GET_VAR( ncid, ncalvid , src_calving (:,:) )
  173. nret = NF90_GET_VAR( ncid, ncalvhid, src_calving_hflx (:,:) )
  174. nret = NF90_GET_VAR( ncid, nsheatid, berg_grid%stored_heat(:,:) )
  175. nret = NF90_GET_VAR( ncid, nkountid, idata2(:) )
  176. num_bergs(:) = idata2(:)
  177. ! Finish up
  178. nret = NF90_CLOSE(ncid)
  179. IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart: nf_close failed')
  180. ! Sanity check
  181. jn = icb_utl_count()
  182. IF (nn_verbose_level >= 0) &
  183. WRITE(numout,'(2(a,i5))') 'icebergs, read_restart_bergs: # bergs =',jn,' on PE',narea-1
  184. IF( lk_mpp ) THEN
  185. CALL mpp_sum(ibergs_in_file)
  186. CALL mpp_sum(jn)
  187. ENDIF
  188. IF(lwp) WRITE(numout,'(a,i5,a,i5,a)') 'icebergs, read_restart_bergs: there were',ibergs_in_file, &
  189. & ' bergs in the restart file and', jn,' bergs have been read'
  190. !
  191. IF( lwp .and. nn_verbose_level >= 0) WRITE(numout,'(a)') 'icebergs, read_restart_bergs: completed'
  192. !
  193. END SUBROUTINE icb_rst_read
  194. SUBROUTINE icb_rst_write( kt )
  195. !!----------------------------------------------------------------------
  196. !! *** SUBROUTINE icb_rst_write ***
  197. !!
  198. !!----------------------------------------------------------------------
  199. INTEGER, INTENT( in ) :: kt
  200. !
  201. INTEGER :: jn ! dummy loop index
  202. INTEGER :: ix_dim, iy_dim, ik_dim, in_dim
  203. CHARACTER(len=256) :: cl_path
  204. CHARACTER(len=256) :: cl_filename
  205. TYPE(iceberg), POINTER :: this
  206. TYPE(point) , POINTER :: pt
  207. !!----------------------------------------------------------------------
  208. ! Assume we write iceberg restarts to same directory as ocean restarts.
  209. cl_path = TRIM(cn_ocerst_outdir)
  210. IF( cl_path(LEN_TRIM(cl_path):) /= '/' ) cl_path = TRIM(cl_path) // '/'
  211. IF( lk_mpp ) THEN
  212. WRITE(cl_filename,'(A,"_icebergs_",I8.8,"_restart_",I4.4,".nc")') TRIM(cexper), kt, narea-1
  213. ELSE
  214. WRITE(cl_filename,'(A,"_icebergs_",I8.8,"_restart.nc")') TRIM(cexper), kt
  215. ENDIF
  216. IF (nn_verbose_level >= 0) WRITE(numout,'(2a)') 'icebergs, write_restart: creating ',TRIM(cl_path)//TRIM(cl_filename)
  217. nret = NF90_CREATE(TRIM(cl_path)//TRIM(cl_filename), NF90_CLOBBER, ncid)
  218. IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_create failed')
  219. ! Dimensions
  220. nret = NF90_DEF_DIM(ncid, 'x', jpi, ix_dim)
  221. IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim x failed')
  222. nret = NF90_DEF_DIM(ncid, 'y', jpj, iy_dim)
  223. IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim y failed')
  224. nret = NF90_DEF_DIM(ncid, 'c', nclasses, nc_dim)
  225. IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim c failed')
  226. nret = NF90_DEF_DIM(ncid, 'k', nkounts, ik_dim)
  227. IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim k failed')
  228. ! global attributes
  229. IF( lk_mpp ) THEN
  230. ! Set domain parameters (assume jpdom_local_full)
  231. nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_number_total' , jpnij )
  232. nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_number' , narea-1 )
  233. nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_dimensions_ids' , (/1 , 2 /) )
  234. nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_size_global' , (/jpiglo, jpjglo/) )
  235. nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_size_local' , (/jpi , jpj /) )
  236. nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_position_first' , (/nimpp , njmpp /) )
  237. nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_position_last' , (/nimpp + jpi - 1 , njmpp + jpj - 1 /) )
  238. nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_halo_size_start', (/nldi - 1 , nldj - 1 /) )
  239. nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_halo_size_end' , (/jpi - nlei , jpj - nlej /) )
  240. nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_type' , 'BOX' )
  241. ENDIF
  242. IF (associated(first_berg)) then
  243. nret = NF90_DEF_DIM(ncid, 'n', NF90_UNLIMITED, in_dim)
  244. IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim n failed')
  245. ENDIF
  246. ! Variables
  247. nret = NF90_DEF_VAR(ncid, 'kount' , NF90_INT , (/ ik_dim /), nkountid)
  248. nret = NF90_DEF_VAR(ncid, 'calving' , NF90_DOUBLE, (/ ix_dim, iy_dim /), ncalvid)
  249. nret = NF90_DEF_VAR(ncid, 'calving_hflx', NF90_DOUBLE, (/ ix_dim, iy_dim /), ncalvhid)
  250. nret = NF90_DEF_VAR(ncid, 'stored_ice' , NF90_DOUBLE, (/ ix_dim, iy_dim, nc_dim /), nsiceid)
  251. nret = NF90_DEF_VAR(ncid, 'stored_heat' , NF90_DOUBLE, (/ ix_dim, iy_dim /), nsheatid)
  252. ! Attributes
  253. nret = NF90_PUT_ATT(ncid, ncalvid , 'long_name', 'iceberg calving')
  254. nret = NF90_PUT_ATT(ncid, ncalvid , 'units', 'some')
  255. nret = NF90_PUT_ATT(ncid, ncalvhid, 'long_name', 'heat flux associated with iceberg calving')
  256. nret = NF90_PUT_ATT(ncid, ncalvhid, 'units', 'some')
  257. nret = NF90_PUT_ATT(ncid, nsiceid , 'long_name', 'stored ice used to calve icebergs')
  258. nret = NF90_PUT_ATT(ncid, nsiceid , 'units', 'kg/s')
  259. nret = NF90_PUT_ATT(ncid, nsheatid, 'long_name', 'heat in stored ice used to calve icebergs')
  260. nret = NF90_PUT_ATT(ncid, nsheatid, 'units', 'J/kg/s')
  261. IF ( ASSOCIATED(first_berg) ) THEN
  262. ! Only add berg variables for this PE if we have anything to say
  263. ! Variables
  264. nret = NF90_DEF_VAR(ncid, 'lon', NF90_DOUBLE, in_dim, nlonid)
  265. nret = NF90_DEF_VAR(ncid, 'lat', NF90_DOUBLE, in_dim, nlatid)
  266. nret = NF90_DEF_VAR(ncid, 'xi', NF90_DOUBLE, in_dim, nxid)
  267. nret = NF90_DEF_VAR(ncid, 'yj', NF90_DOUBLE, in_dim, nyid)
  268. nret = NF90_DEF_VAR(ncid, 'uvel', NF90_DOUBLE, in_dim, nuvelid)
  269. nret = NF90_DEF_VAR(ncid, 'vvel', NF90_DOUBLE, in_dim, nvvelid)
  270. nret = NF90_DEF_VAR(ncid, 'mass', NF90_DOUBLE, in_dim, nmassid)
  271. nret = NF90_DEF_VAR(ncid, 'thickness', NF90_DOUBLE, in_dim, nthicknessid)
  272. nret = NF90_DEF_VAR(ncid, 'width', NF90_DOUBLE, in_dim, nwidthid)
  273. nret = NF90_DEF_VAR(ncid, 'length', NF90_DOUBLE, in_dim, nlengthid)
  274. nret = NF90_DEF_VAR(ncid, 'number', NF90_INT, (/ik_dim,in_dim/), numberid)
  275. nret = NF90_DEF_VAR(ncid, 'year', NF90_INT, in_dim, nyearid)
  276. nret = NF90_DEF_VAR(ncid, 'day', NF90_DOUBLE, in_dim, ndayid)
  277. nret = NF90_DEF_VAR(ncid, 'mass_scaling', NF90_DOUBLE, in_dim, nscaling_id)
  278. nret = NF90_DEF_VAR(ncid, 'mass_of_bits', NF90_DOUBLE, in_dim, nmass_of_bits_id)
  279. nret = NF90_DEF_VAR(ncid, 'heat_density', NF90_DOUBLE, in_dim, nheat_density_id)
  280. ! Attributes
  281. nret = NF90_PUT_ATT(ncid, nlonid, 'long_name', 'longitude')
  282. nret = NF90_PUT_ATT(ncid, nlonid, 'units', 'degrees_E')
  283. nret = NF90_PUT_ATT(ncid, nlatid, 'long_name', 'latitude')
  284. nret = NF90_PUT_ATT(ncid, nlatid, 'units', 'degrees_N')
  285. nret = NF90_PUT_ATT(ncid, nxid, 'long_name', 'x grid box position')
  286. nret = NF90_PUT_ATT(ncid, nxid, 'units', 'fractional')
  287. nret = NF90_PUT_ATT(ncid, nyid, 'long_name', 'y grid box position')
  288. nret = NF90_PUT_ATT(ncid, nyid, 'units', 'fractional')
  289. nret = NF90_PUT_ATT(ncid, nuvelid, 'long_name', 'zonal velocity')
  290. nret = NF90_PUT_ATT(ncid, nuvelid, 'units', 'm/s')
  291. nret = NF90_PUT_ATT(ncid, nvvelid, 'long_name', 'meridional velocity')
  292. nret = NF90_PUT_ATT(ncid, nvvelid, 'units', 'm/s')
  293. nret = NF90_PUT_ATT(ncid, nmassid, 'long_name', 'mass')
  294. nret = NF90_PUT_ATT(ncid, nmassid, 'units', 'kg')
  295. nret = NF90_PUT_ATT(ncid, nthicknessid, 'long_name', 'thickness')
  296. nret = NF90_PUT_ATT(ncid, nthicknessid, 'units', 'm')
  297. nret = NF90_PUT_ATT(ncid, nwidthid, 'long_name', 'width')
  298. nret = NF90_PUT_ATT(ncid, nwidthid, 'units', 'm')
  299. nret = NF90_PUT_ATT(ncid, nlengthid, 'long_name', 'length')
  300. nret = NF90_PUT_ATT(ncid, nlengthid, 'units', 'm')
  301. nret = NF90_PUT_ATT(ncid, numberid, 'long_name', 'iceberg number on this processor')
  302. nret = NF90_PUT_ATT(ncid, numberid, 'units', 'count')
  303. nret = NF90_PUT_ATT(ncid, nyearid, 'long_name', 'calendar year of calving event')
  304. nret = NF90_PUT_ATT(ncid, nyearid, 'units', 'years')
  305. nret = NF90_PUT_ATT(ncid, ndayid, 'long_name', 'year day of calving event')
  306. nret = NF90_PUT_ATT(ncid, ndayid, 'units', 'days')
  307. nret = NF90_PUT_ATT(ncid, nscaling_id, 'long_name', 'scaling factor for mass of calving berg')
  308. nret = NF90_PUT_ATT(ncid, nscaling_id, 'units', 'none')
  309. nret = NF90_PUT_ATT(ncid, nmass_of_bits_id, 'long_name', 'mass of bergy bits')
  310. nret = NF90_PUT_ATT(ncid, nmass_of_bits_id, 'units', 'kg')
  311. nret = NF90_PUT_ATT(ncid, nheat_density_id, 'long_name', 'heat density')
  312. nret = NF90_PUT_ATT(ncid, nheat_density_id, 'units', 'J/kg')
  313. ENDIF ! associated(first_berg)
  314. ! End define mode
  315. nret = NF90_ENDDEF(ncid)
  316. ! --------------------------------
  317. ! now write some data
  318. nstrt3(1) = 1
  319. nstrt3(2) = 1
  320. nlngth3(1) = jpi
  321. nlngth3(2) = jpj
  322. nlngth3(3) = 1
  323. DO jn=1,nclasses
  324. griddata(:,:,1) = berg_grid%stored_ice(:,:,jn)
  325. nstrt3(3) = jn
  326. nret = NF90_PUT_VAR( ncid, nsiceid, griddata, nstrt3, nlngth3 )
  327. IF (nret .ne. NF90_NOERR) THEN
  328. IF( lwp ) WRITE(numout,*) TRIM(NF90_STRERROR( nret ))
  329. CALL ctl_stop('icebergs, write_restart: nf_put_var stored_ice failed')
  330. ENDIF
  331. ENDDO
  332. IF( lwp ) WRITE(numout,*) 'file: ',TRIM(cl_path)//TRIM(cl_filename),' var: stored_ice written'
  333. nret = NF90_PUT_VAR( ncid, nkountid, num_bergs(:) )
  334. IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var kount failed')
  335. nret = NF90_PUT_VAR( ncid, nsheatid, berg_grid%stored_heat(:,:) )
  336. IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var stored_heat failed')
  337. IF( lwp ) WRITE(numout,*) 'file: ',TRIM(cl_path)//TRIM(cl_filename),' var: stored_heat written'
  338. nret = NF90_PUT_VAR( ncid, ncalvid , src_calving(:,:) )
  339. IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var calving failed')
  340. nret = NF90_PUT_VAR( ncid, ncalvhid, src_calving_hflx(:,:) )
  341. IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var calving_hflx failed')
  342. IF( lwp ) WRITE(numout,*) 'file: ',TRIM(cl_path)//TRIM(cl_filename),' var: calving written'
  343. IF ( ASSOCIATED(first_berg) ) THEN
  344. ! Write variables
  345. ! just write out the current point of the trajectory
  346. this => first_berg
  347. jn = 0
  348. DO WHILE (ASSOCIATED(this))
  349. pt => this%current_point
  350. jn=jn+1
  351. nret = NF90_PUT_VAR(ncid, numberid, this%number, (/1,jn/), (/nkounts,1/) )
  352. nret = NF90_PUT_VAR(ncid, nscaling_id, this%mass_scaling, (/ jn /) )
  353. nret = NF90_PUT_VAR(ncid, nlonid, pt%lon, (/ jn /) )
  354. nret = NF90_PUT_VAR(ncid, nlatid, pt%lat, (/ jn /) )
  355. nret = NF90_PUT_VAR(ncid, nxid, pt%xi, (/ jn /) )
  356. nret = NF90_PUT_VAR(ncid, nyid, pt%yj, (/ jn /) )
  357. nret = NF90_PUT_VAR(ncid, nuvelid, pt%uvel, (/ jn /) )
  358. nret = NF90_PUT_VAR(ncid, nvvelid, pt%vvel, (/ jn /) )
  359. nret = NF90_PUT_VAR(ncid, nmassid, pt%mass, (/ jn /) )
  360. nret = NF90_PUT_VAR(ncid, nthicknessid, pt%thickness, (/ jn /) )
  361. nret = NF90_PUT_VAR(ncid, nwidthid, pt%width, (/ jn /) )
  362. nret = NF90_PUT_VAR(ncid, nlengthid, pt%length, (/ jn /) )
  363. nret = NF90_PUT_VAR(ncid, nyearid, pt%year, (/ jn /) )
  364. nret = NF90_PUT_VAR(ncid, ndayid, pt%day, (/ jn /) )
  365. nret = NF90_PUT_VAR(ncid, nmass_of_bits_id, pt%mass_of_bits, (/ jn /) )
  366. nret = NF90_PUT_VAR(ncid, nheat_density_id, pt%heat_density, (/ jn /) )
  367. this=>this%next
  368. END DO
  369. !
  370. ENDIF ! associated(first_berg)
  371. ! Finish up
  372. nret = NF90_CLOSE(ncid)
  373. IF (nret /= NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_close failed')
  374. !
  375. END SUBROUTINE icb_rst_write
  376. !
  377. END MODULE icbrst