Changeset 14078 for NEMO/branches/UKMO/NEMO_4.0.4_GO8_package/src
- Timestamp:
- 2020-12-04T11:17:01+01:00 (4 years ago)
- Location:
- NEMO/branches/UKMO/NEMO_4.0.4_GO8_package/src
- Files:
-
- 14 edited
- 1 copied
Legend:
- Unmodified
- Added
- Removed
-
NEMO/branches/UKMO/NEMO_4.0.4_GO8_package/src/ICE/icerst.F90
r14075 r14078 27 27 USE in_out_manager ! I/O manager 28 28 USE iom ! I/O manager library 29 USE ioipsl, ONLY : ju2ymds ! for calendar 29 30 USE lib_mpp ! MPP library 30 31 USE lib_fortran ! fortran utilities (glob_sum + no signed zero) … … 52 53 INTEGER, INTENT(in) :: kt ! number of iteration 53 54 ! 55 INTEGER :: iyear, imonth, iday 56 REAL (wp) :: zsec 57 REAL (wp) :: zfjulday 54 58 CHARACTER(len=20) :: clkt ! ocean time-step define as a character 55 59 CHARACTER(len=50) :: clname ! ice output restart file name … … 67 71 IF( nitrst <= nitend .AND. nitrst > 0 ) THEN 68 72 ! beware of the format used to write kt (default is i8.8, that should be large enough...) 69 IF( nitrst > 99999999 ) THEN ; WRITE(clkt, * ) nitrst 70 ELSE ; WRITE(clkt, '(i8.8)') nitrst 73 IF ( ln_rstdate ) THEN 74 zfjulday = fjulday + (2*nn_fsbc+1)*rdt / rday 75 IF( ABS(zfjulday - REAL(NINT(zfjulday),wp)) < 0.1 / rday ) zfjulday = REAL(NINT(zfjulday),wp) ! avoid truncation error 76 CALL ju2ymds( zfjulday, iyear, imonth, iday, zsec ) 77 WRITE(clkt, '(i4.4,2i2.2)') iyear, imonth, iday 78 ELSE 79 IF( nitrst > 99999999 ) THEN ; WRITE(clkt, * ) nitrst 80 ELSE ; WRITE(clkt, '(i8.8)') nitrst 81 ENDIF 71 82 ENDIF 72 83 ! create the file -
NEMO/branches/UKMO/NEMO_4.0.4_GO8_package/src/ICE/icewri.F90
r14075 r14078 104 104 IF( iom_use('snwthic' ) ) CALL iom_put( 'snwthic', hm_s * zmsk00 ) ! snw thickness 105 105 IF( iom_use('icebrv' ) ) CALL iom_put( 'icebrv' , bvm_i* 100. * zmsk00 ) ! brine volume 106 IF( iom_use('iceage' ) ) CALL iom_put( 'iceage' , om_i / rday * zmsk15 + zmiss_val * ( 1._wp - zmsk15 ) )! ice age106 IF( iom_use('iceage' ) ) CALL iom_put( 'iceage' , om_i / rday * zmsk15 ) ! ice age 107 107 IF( iom_use('icehnew' ) ) CALL iom_put( 'icehnew', ht_i_new ) ! new ice thickness formed in the leads 108 108 IF( iom_use('snwvolu' ) ) CALL iom_put( 'snwvolu', vt_s * zmsksn ) ! snow volume … … 119 119 IF( iom_use('icevlid' ) ) CALL iom_put( 'icevlid', vt_il * zmsk00 ) ! melt pond lid total volume per unit area 120 120 ! salt 121 IF( iom_use('icesalt' ) ) CALL iom_put( 'icesalt', sm_i * zmsk00 + zmiss_val * ( 1._wp - zmsk00 ) )! mean ice salinity121 IF( iom_use('icesalt' ) ) CALL iom_put( 'icesalt', sm_i * zmsk00 ) ! mean ice salinity 122 122 IF( iom_use('icesalm' ) ) CALL iom_put( 'icesalm', st_i * rhoi * 1.0e-3 * zmsk00 ) ! Mass of salt in sea ice per cell area 123 123 ! heat 124 IF( iom_use('icetemp' ) ) CALL iom_put( 'icetemp', ( tm_i - rt0 ) * zmsk00 + zmiss_val * ( 1._wp - zmsk00 ) )! ice mean temperature125 IF( iom_use('snwtemp' ) ) CALL iom_put( 'snwtemp', ( tm_s - rt0 ) * zmsksn + zmiss_val * ( 1._wp - zmsksn ) )! snw mean temperature126 IF( iom_use('icettop' ) ) CALL iom_put( 'icettop', ( tm_su - rt0 ) * zmsk00 + zmiss_val * ( 1._wp - zmsk00 ) )! temperature at the ice surface127 IF( iom_use('icetbot' ) ) CALL iom_put( 'icetbot', ( t_bo - rt0 ) * zmsk00 + zmiss_val * ( 1._wp - zmsk00 ) )! temperature at the ice bottom128 IF( iom_use('icetsni' ) ) CALL iom_put( 'icetsni', ( tm_si - rt0 ) * zmsk00 + zmiss_val * ( 1._wp - zmsk00 ) )! temperature at the snow-ice interface124 IF( iom_use('icetemp' ) ) CALL iom_put( 'icetemp', ( tm_i - rt0 ) * zmsk00 ) ! ice mean temperature 125 IF( iom_use('snwtemp' ) ) CALL iom_put( 'snwtemp', ( tm_s - rt0 ) * zmsksn ) ! snw mean temperature 126 IF( iom_use('icettop' ) ) CALL iom_put( 'icettop', ( tm_su - rt0 ) * zmsk00 ) ! temperature at the ice surface 127 IF( iom_use('icetbot' ) ) CALL iom_put( 'icetbot', ( t_bo - rt0 ) * zmsk00 ) ! temperature at the ice bottom 128 IF( iom_use('icetsni' ) ) CALL iom_put( 'icetsni', ( tm_si - rt0 ) * zmsk00 ) ! temperature at the snow-ice interface 129 129 IF( iom_use('icehc' ) ) CALL iom_put( 'icehc' , -et_i * zmsk00 ) ! ice heat content 130 130 IF( iom_use('snwhc' ) ) CALL iom_put( 'snwhc' , -et_s * zmsksn ) ! snow heat content … … 153 153 IF( iom_use('icemask_cat' ) ) CALL iom_put( 'icemask_cat' , zmsk00l ) ! ice mask 0% 154 154 IF( iom_use('iceconc_cat' ) ) CALL iom_put( 'iceconc_cat' , a_i * zmsk00l ) ! area for categories 155 IF( iom_use('icethic_cat' ) ) CALL iom_put( 'icethic_cat' , h_i * zmsk00l + zmiss_val * ( 1._wp - zmsk00l )) ! thickness for categories156 IF( iom_use('snwthic_cat' ) ) CALL iom_put( 'snwthic_cat' , h_s * zmsksnl + zmiss_val * ( 1._wp - zmsksnl )) ! snow depth for categories157 IF( iom_use('icesalt_cat' ) ) CALL iom_put( 'icesalt_cat' , s_i * zmsk00l + zmiss_val * ( 1._wp - zmsk00l )) ! salinity for categories158 IF( iom_use('iceage_cat' ) ) CALL iom_put( 'iceage_cat' , o_i / rday * zmsk00l + zmiss_val * ( 1._wp - zmsk00l )) ! ice age155 IF( iom_use('icethic_cat' ) ) CALL iom_put( 'icethic_cat' , h_i * zmsk00l ) ! thickness for categories 156 IF( iom_use('snwthic_cat' ) ) CALL iom_put( 'snwthic_cat' , h_s * zmsksnl ) ! snow depth for categories 157 IF( iom_use('icesalt_cat' ) ) CALL iom_put( 'icesalt_cat' , s_i * zmsk00l ) ! salinity for categories 158 IF( iom_use('iceage_cat' ) ) CALL iom_put( 'iceage_cat' , o_i / rday * zmsk00l ) ! ice age 159 159 IF( iom_use('icetemp_cat' ) ) CALL iom_put( 'icetemp_cat' , ( SUM( t_i, dim=3 ) * r1_nlay_i - rt0 ) & 160 & * zmsk00l + zmiss_val * ( 1._wp - zmsk00l )) ! ice temperature160 & * zmsk00l ) ! ice temperature 161 161 IF( iom_use('snwtemp_cat' ) ) CALL iom_put( 'snwtemp_cat' , ( SUM( t_s, dim=3 ) * r1_nlay_s - rt0 ) & 162 & * zmsksnl + zmiss_val * ( 1._wp - zmsksnl )) ! snow temperature163 IF( iom_use('icettop_cat' ) ) CALL iom_put( 'icettop_cat' , ( t_su - rt0 ) * zmsk00l + zmiss_val * ( 1._wp - zmsk00l )) ! surface temperature164 IF( iom_use('icebrv_cat' ) ) CALL iom_put( 'icebrv_cat' , bv_i * 100. * zmsk00l + zmiss_val * ( 1._wp - zmsk00l )) ! brine volume162 & * zmsksnl ) ! snow temperature 163 IF( iom_use('icettop_cat' ) ) CALL iom_put( 'icettop_cat' , ( t_su - rt0 ) * zmsk00l ) ! surface temperature 164 IF( iom_use('icebrv_cat' ) ) CALL iom_put( 'icebrv_cat' , bv_i * 100. * zmsk00l ) ! brine volume 165 165 IF( iom_use('iceapnd_cat' ) ) CALL iom_put( 'iceapnd_cat' , a_ip * zmsk00l ) ! melt pond frac for categories 166 166 IF( iom_use('icehpnd_cat' ) ) CALL iom_put( 'icehpnd_cat' , h_ip * zmsk00l + zmiss_val * ( 1._wp - zmsk00l ) ) ! melt pond thickness for categories -
NEMO/branches/UKMO/NEMO_4.0.4_GO8_package/src/OCE/DIA/diawri.F90
r14075 r14078 49 49 USE iom ! 50 50 USE ioipsl ! 51 51 USE eosbn2 52 52 #if defined key_si3 53 53 USE ice … … 113 113 REAL(wp), DIMENSION(jpi,jpj) :: z2d ! 2D workspace 114 114 REAL(wp), DIMENSION(jpi,jpj,jpk) :: z3d ! 3D workspace 115 CHARACTER(len=4),SAVE :: ttype , stype ! temperature and salinity type 115 116 !!---------------------------------------------------------------------- 116 117 ! 118 IF( kt == nit000 ) THEN 119 IF( ln_TEOS10 ) THEN 120 IF ( iom_use("toce_pot") .OR. iom_use("soce_pra") .OR. iom_use("sst_pot") .OR. iom_use("sss_pra") & 121 & .OR. iom_use("sbt_pot") .OR. iom_use("sbs_pra") .OR. iom_use("sstgrad_pot") .OR. iom_use("sstgrad2_pot") & 122 & .OR. iom_use("tosmint_pot") .OR. iom_use("somint_pra")) THEN 123 CALL ctl_stop( 'diawri: potential temperature and practical salinity not available with ln_TEOS10' ) 124 ELSE 125 ttype='con' ; stype='abs' ! teos-10 using conservative temperature and absolute salinity 126 ENDIF 127 ELSE IF( ln_EOS80 ) THEN 128 IF ( iom_use("toce_con") .OR. iom_use("soce_abs") .OR. iom_use("sst_con") .OR. iom_use("sss_abs") & 129 & .OR. iom_use("sbt_con") .OR. iom_use("sbs_abs") .OR. iom_use("sstgrad_con") .OR. iom_use("sstgrad2_con") & 130 & .OR. iom_use("tosmint_con") .OR. iom_use("somint_abs")) THEN 131 CALL ctl_stop( 'diawri: conservative temperature and absolute salinity not available with ln_EOS80' ) 132 ELSE 133 ttype='pot' ; stype='pra' ! eos-80 using potential temperature and practical salinity 134 ENDIF 135 ELSE IF ( ln_SEOS) THEN 136 ttype='seos' ; stype='seos' ! seos using Simplified Equation of state 137 ENDIF 138 ENDIF 139 117 140 IF( ln_timing ) CALL timing_start('dia_wri') 118 141 ! … … 144 167 CALL iom_put( "wetdep" , ht_0(:,:) + sshn(:,:) ) 145 168 146 CALL iom_put( "toce ", tsn(:,:,:,jp_tem) ) ! 3D temperature147 CALL iom_put( "sst ", tsn(:,:,1,jp_tem) ) ! surface temperature148 IF ( iom_use("sbt ") ) THEN169 CALL iom_put( "toce_"//ttype, tsn(:,:,:,jp_tem) ) ! 3D temperature 170 CALL iom_put( "sst_"//ttype, tsn(:,:,1,jp_tem) ) ! surface temperature 171 IF ( iom_use("sbt_"//ttype) ) THEN 149 172 DO jj = 1, jpj 150 173 DO ji = 1, jpi … … 153 176 END DO 154 177 END DO 155 CALL iom_put( "sbt ", z2d ) ! bottom temperature178 CALL iom_put( "sbt_"//ttype, z2d ) ! bottom temperature 156 179 ENDIF 157 180 158 CALL iom_put( "soce ", tsn(:,:,:,jp_sal) ) ! 3D salinity159 CALL iom_put( "sss ", tsn(:,:,1,jp_sal) ) ! surface salinity160 IF ( iom_use("sbs ") ) THEN181 CALL iom_put( "soce_"//stype, tsn(:,:,:,jp_sal) ) ! 3D salinity 182 CALL iom_put( "sss_"//stype, tsn(:,:,1,jp_sal) ) ! surface salinity 183 IF ( iom_use("sbs_"//stype) ) THEN 161 184 DO jj = 1, jpj 162 185 DO ji = 1, jpi … … 165 188 END DO 166 189 END DO 167 CALL iom_put( "sbs ", z2d ) ! bottom salinity190 CALL iom_put( "sbs_"//stype, z2d ) ! bottom salinity 168 191 ENDIF 169 192 … … 233 256 IF( iom_use('logavs') ) CALL iom_put( "logavs", LOG( MAX( 1.e-20_wp, avs(:,:,:) ) ) ) 234 257 235 IF ( iom_use("sstgrad ") .OR. iom_use("sstgrad2") ) THEN258 IF ( iom_use("sstgrad_"//ttype) .OR. iom_use("sstgrad2_"//ttype) ) THEN 236 259 DO jj = 2, jpjm1 ! sst gradient 237 260 DO ji = fs_2, fs_jpim1 ! vector opt. … … 244 267 END DO 245 268 CALL lbc_lnk( 'diawri', z2d, 'T', 1. ) 246 CALL iom_put( "sstgrad2 ", z2d ) ! square of module of sst gradient269 CALL iom_put( "sstgrad2_"//ttype, z2d ) ! square of module of sst gradient 247 270 z2d(:,:) = SQRT( z2d(:,:) ) 248 CALL iom_put( "sstgrad ", z2d ) ! module of sst gradient271 CALL iom_put( "sstgrad_"//ttype , z2d ) ! module of sst gradient 249 272 ENDIF 250 273 … … 365 388 ENDIF 366 389 367 IF( iom_use("tosmint ") ) THEN390 IF( iom_use("tosmint_"//ttype) ) THEN 368 391 z2d(:,:) = 0._wp 369 392 DO jk = 1, jpkm1 … … 375 398 END DO 376 399 CALL lbc_lnk( 'diawri', z2d, 'T', -1. ) 377 CALL iom_put( "tosmint ", rau0 * z2d ) ! Vertical integral of temperature378 ENDIF 379 IF( iom_use("somint ") ) THEN400 CALL iom_put( "tosmint_"//ttype, rau0 * z2d ) ! Vertical integral of temperature 401 ENDIF 402 IF( iom_use("somint_"//stype) ) THEN 380 403 z2d(:,:)=0._wp 381 404 DO jk = 1, jpkm1 … … 387 410 END DO 388 411 CALL lbc_lnk( 'diawri', z2d, 'T', -1. ) 389 CALL iom_put( "somint ", rau0 * z2d ) ! Vertical integral of salinity412 CALL iom_put( "somint_"//stype, rau0 * z2d ) ! Vertical integral of salinity 390 413 ENDIF 391 414 … … 930 953 CALL iom_rstput( 0, 0, inum, 'sdvecrtz', wsd ) ! now StokesDrift k-velocity 931 954 ENDIF 932 933 955 #if defined key_si3 934 956 IF( nn_ice == 2 ) THEN ! condition needed in case agrif + ice-model but no-ice in child grid -
NEMO/branches/UKMO/NEMO_4.0.4_GO8_package/src/OCE/DOM/domain.F90
r14075 r14078 292 292 & nn_it000, nn_itend , nn_date0 , nn_time0 , nn_leapy , nn_istate , & 293 293 & nn_stock, nn_write , ln_mskland , ln_clobber , nn_chunksz, nn_euler , & 294 & ln_cfmeta, ln_iscpl, ln_xios_read, nn_wxios 294 & ln_cfmeta, ln_iscpl, ln_xios_read, nn_wxios, ln_rstdate, ln_rst_eos 295 295 296 NAMELIST/namdom/ ln_linssh, rn_isfhmin, rn_rdt, rn_atfp, ln_crs, ln_meshmask 296 297 #if defined key_netcdf4 … … 340 341 #endif 341 342 WRITE(numout,*) ' mask land points ln_mskland = ', ln_mskland 343 WRITE(numout,*) ' date-stamp restart files ln_rstdate = ', ln_rstdate 342 344 WRITE(numout,*) ' additional CF standard metadata ln_cfmeta = ', ln_cfmeta 343 345 WRITE(numout,*) ' overwrite an existing file ln_clobber = ', ln_clobber 344 346 WRITE(numout,*) ' NetCDF chunksize (bytes) nn_chunksz = ', nn_chunksz 345 347 WRITE(numout,*) ' IS coupling at the restart step ln_iscpl = ', ln_iscpl 348 WRITE(numout,*) ' check restart equation of state ln_rst_eos = ', ln_rst_eos 349 346 350 IF( TRIM(Agrif_CFixed()) == '0' ) THEN 347 351 WRITE(numout,*) ' READ restart for a single file using XIOS ln_xios_read =', ln_xios_read -
NEMO/branches/UKMO/NEMO_4.0.4_GO8_package/src/OCE/DOM/dommsk.F90
r14075 r14078 32 32 USE lbclnk ! ocean lateral boundary conditions (or mpp link) 33 33 USE lib_mpp ! Massively Parallel Processing library 34 USE iom ! For shlat2d 35 USE fldread ! for sn_shlat2d 34 36 35 37 IMPLICIT NONE … … 93 95 INTEGER :: ios, inum 94 96 !! 95 NAMELIST/namlbc/ rn_shlat, ln_vorlat 97 REAL(wp) :: zshlat !: locally modified shlat for some strait 98 REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zshlat2d 99 LOGICAL :: ln_shlat2d 100 CHARACTER(len = 256) :: cn_shlat2d_file, cn_shlat2d_var 101 !! 102 NAMELIST/namlbc/ rn_shlat, ln_vorlat, ln_shlat2d, cn_shlat2d_file, cn_shlat2d_var 96 103 NAMELIST/nambdy/ ln_bdy ,nb_bdy, ln_coords_file, cn_coords_file, & 97 104 & ln_mask_file, cn_mask_file, cn_dyn2d, nn_dyn2d_dta, & … … 120 127 ! 121 128 IF(lwp) WRITE(numout,*) 122 IF ( rn_shlat == 0. ) THEN ; IF(lwp) WRITE(numout,*) ' ==>>> ocean lateral free-slip' 123 ELSEIF ( rn_shlat == 2. ) THEN ; IF(lwp) WRITE(numout,*) ' ==>>> ocean lateral no-slip' 124 ELSEIF ( 0. < rn_shlat .AND. rn_shlat < 2. ) THEN ; IF(lwp) WRITE(numout,*) ' ==>>> ocean lateral partial-slip' 125 ELSEIF ( 2. < rn_shlat ) THEN ; IF(lwp) WRITE(numout,*) ' ==>>> ocean lateral strong-slip' 129 130 IF ( ln_shlat2d ) THEN 131 IF(lwp) WRITE(numout,*) ' READ shlat as a 2D coefficient in a file ' 132 ALLOCATE( zshlat2d(jpi,jpj) ) 133 CALL iom_open(TRIM(cn_shlat2d_file), inum) 134 CALL iom_get (inum, jpdom_data, TRIM(cn_shlat2d_var), zshlat2d, 1) ! 135 CALL iom_close(inum) 126 136 ELSE 127 CALL ctl_stop( 'dom_msk: wrong value for rn_shlat (i.e. a negalive value). We stop.' ) 137 IF ( rn_shlat == 0. ) THEN ; IF(lwp) WRITE(numout,*) ' ==>>> ocean lateral free-slip' 138 ELSEIF ( rn_shlat == 2. ) THEN ; IF(lwp) WRITE(numout,*) ' ==>>> ocean lateral no-slip' 139 ELSEIF ( 0. < rn_shlat .AND. rn_shlat < 2. ) THEN ; IF(lwp) WRITE(numout,*) ' ==>>> ocean lateral partial-slip' 140 ELSEIF ( 2. < rn_shlat ) THEN ; IF(lwp) WRITE(numout,*) ' ==>>> ocean lateral strong-slip' 141 ELSE 142 CALL ctl_stop( 'dom_msk: wrong value for rn_shlat (i.e. a negalive value). We stop.' ) 143 ENDIF 128 144 ENDIF 129 145 … … 240 256 ! Lateral boundary conditions on velocity (modify fmask) 241 257 ! --------------------------------------- 242 IF( rn_shlat /= 0 ) THEN ! Not free-slip lateral boundary condition258 IF( rn_shlat /= 0 .or. ln_shlat2d ) THEN ! Not free-slip lateral boundary condition everywhere 243 259 ! 244 260 DO jk = 1, jpk 245 DO jj = 2, jpjm1 246 DO ji = fs_2, fs_jpim1 ! vector opt. 247 IF( fmask(ji,jj,jk) == 0._wp ) THEN 248 fmask(ji,jj,jk) = rn_shlat * MIN( 1._wp , MAX( umask(ji,jj,jk), umask(ji,jj+1,jk), & 249 & vmask(ji,jj,jk), vmask(ji+1,jj,jk) ) ) 250 ENDIF 261 IF ( ln_shlat2d ) THEN 262 DO jj = 2, jpjm1 263 DO ji = fs_2, fs_jpim1 ! vector opt. 264 IF( fmask(ji,jj,jk) == 0._wp ) THEN 265 fmask(ji,jj,jk) = zshlat2d(ji,jj) * MIN( 1._wp , MAX( umask(ji,jj,jk), umask(ji,jj+1,jk), & 266 & vmask(ji,jj,jk), vmask(ji+1,jj,jk) ) ) 267 ENDIF 268 END DO 251 269 END DO 252 END DO 270 ELSE 271 DO jj = 2, jpjm1 272 DO ji = fs_2, fs_jpim1 ! vector opt. 273 IF( fmask(ji,jj,jk) == 0._wp ) THEN 274 fmask(ji,jj,jk) = rn_shlat * MIN( 1._wp , MAX( umask(ji,jj,jk), umask(ji,jj+1,jk), & 275 & vmask(ji,jj,jk), vmask(ji+1,jj,jk) ) ) 276 ENDIF 277 END DO 278 END DO 279 ENDIF 253 280 DO jj = 2, jpjm1 254 281 IF( fmask(1,jj,jk) == 0._wp ) THEN … … 277 304 END DO 278 305 ! 306 IF( ln_shlat2d ) DEALLOCATE( zshlat2d ) 307 ! 279 308 CALL lbc_lnk( 'dommsk', fmask, 'F', 1._wp ) ! Lateral boundary conditions on fmask 280 309 ! … … 284 313 285 314 ! User defined alteration of fmask (use to reduce ocean transport in specified straits) 315 ! Only call if we are not using the shlat2d option. 286 316 ! -------------------------------- 287 317 ! 288 CALL usr_def_fmask( cn_cfg, nn_cfg, fmask ) 318 IF ( .not. ln_shlat2d ) THEN 319 CALL usr_def_fmask( cn_cfg, nn_cfg, fmask ) 320 ENDIF 289 321 ! 290 322 END SUBROUTINE dom_msk -
NEMO/branches/UKMO/NEMO_4.0.4_GO8_package/src/OCE/ICB/icbrst.F90
r14075 r14078 25 25 USE netcdf ! netcdf routines for IO 26 26 USE iom 27 USE ioipsl, ONLY : ju2ymds ! for calendar 27 28 USE icb_oce ! define iceberg arrays 28 29 USE icbutl ! iceberg utility routines … … 191 192 INTEGER :: idg ! number of digits 192 193 INTEGER :: ix_dim, iy_dim, ik_dim, in_dim 193 CHARACTER(len=256) :: cl_path 194 CHARACTER(len=256) :: cl_filename 195 CHARACTER(len=8 ) :: cl_kt 194 INTEGER :: iyear, imonth, iday 195 REAL (wp) :: zsec 196 REAL (wp) :: zfjulday 197 CHARACTER(len=256) :: cl_path 198 CHARACTER(len=256) :: cl_filename 199 CHARACTER(LEN=20) :: clkt ! ocean time-step deine as a character 196 200 CHARACTER(LEN=12 ) :: clfmt ! writing format 197 201 TYPE(iceberg), POINTER :: this … … 209 213 cl_path = TRIM(cn_ocerst_outdir) 210 214 IF( cl_path(LEN_TRIM(cl_path):) /= '/' ) cl_path = TRIM(cl_path) // '/' 211 WRITE(cl_kt, '(i8.8)') kt 212 cl_filename = TRIM(cexper)//"_icebergs_"//cl_kt//"_restart" 215 IF ( ln_rstdate ) THEN 216 zfjulday = fjulday + rdt / rday 217 IF( ABS(zfjulday - REAL(NINT(zfjulday),wp)) < 0.1 / rday ) zfjulday = REAL(NINT(zfjulday),wp) ! avoid truncation error 218 CALL ju2ymds( zfjulday, iyear, imonth, iday, zsec ) 219 WRITE(clkt, '(i4.4,2i2.2)') iyear, imonth, iday 220 ELSE 221 IF( kt > 999999999 ) THEN ; WRITE(clkt, * ) kt 222 ELSE ; WRITE(clkt, '(i8.8)') kt 223 ENDIF 224 ENDIF 225 cl_filename = TRIM(cexper)//"_icebergs_"//TRIM(ADJUSTL(clkt))//"_restart" 213 226 IF( lk_mpp ) THEN 214 227 idg = MAX( INT(LOG10(REAL(MAX(1,jpnij-1),wp))) + 1, 4 ) ! how many digits to we need to write? min=4, max=9 -
NEMO/branches/UKMO/NEMO_4.0.4_GO8_package/src/OCE/IOM/in_out_manager.F90
r14075 r14078 28 28 LOGICAL :: ln_rstart !: start from (F) rest or (T) a restart file 29 29 LOGICAL :: ln_rst_list !: output restarts at list of times (T) or by frequency (F) 30 LOGICAL :: ln_rst_eos !: check equation of state used for the restart is consistent with model 30 31 INTEGER :: nn_rstctl !: control of the time step (0, 1 or 2) 31 32 INTEGER :: nn_rstssh = 0 !: hand made initilization of ssh or not (1/0) … … 40 41 INTEGER, DIMENSION(10) :: nn_stocklist !: restart dump times 41 42 LOGICAL :: ln_mskland !: mask land points in NetCDF outputs (costly: + ~15%) 43 LOGICAL :: ln_rstdate !: T=> stamp output restart files with date instead of timestep 42 44 LOGICAL :: ln_cfmeta !: output additional data to netCDF files required for compliance with the CF metadata standard 43 45 LOGICAL :: ln_clobber !: clobber (overwrite) an existing file -
NEMO/branches/UKMO/NEMO_4.0.4_GO8_package/src/OCE/IOM/iom.F90
r14075 r14078 368 368 !from restart.F90 369 369 CALL iom_set_rstw_var_active("rdt") 370 CALL iom_set_rstw_var_active("neos") 371 370 372 IF ( .NOT. ln_diurnal_only ) THEN 371 373 CALL iom_set_rstw_var_active('ub' ) … … 417 419 i = 0 418 420 i = i + 1; fields(i)%vname="rdt"; fields(i)%grid="grid_scalar" 421 i = i + 1; fields(i)%vname="neos"; fields(i)%grid="grid_scalar" 419 422 i = i + 1; fields(i)%vname="un"; fields(i)%grid="grid_N_3D" 420 423 i = i + 1; fields(i)%vname="ub"; fields(i)%grid="grid_N_3D" -
NEMO/branches/UKMO/NEMO_4.0.4_GO8_package/src/OCE/IOM/restart.F90
r14075 r14078 27 27 USE in_out_manager ! I/O manager 28 28 USE iom ! I/O module 29 USE ioipsl, ONLY : ju2ymds ! for calendar 29 30 USE diurnal_bulk 30 31 USE lib_mpp ! distribued memory computing library … … 59 60 INTEGER, INTENT(in) :: kt ! ocean time-step 60 61 !! 62 INTEGER :: iyear, imonth, iday 63 REAL (wp) :: zsec 64 REAL (wp) :: zfjulday 61 65 CHARACTER(LEN=20) :: clkt ! ocean time-step deine as a character 62 66 CHARACTER(LEN=50) :: clname ! ocean output restart file name … … 90 94 IF( nitrst <= nitend .AND. nitrst > 0 ) THEN 91 95 ! beware of the format used to write kt (default is i8.8, that should be large enough...) 92 IF( nitrst > 999999999 ) THEN ; WRITE(clkt, * ) nitrst 93 ELSE ; WRITE(clkt, '(i8.8)') nitrst 96 IF ( ln_rstdate ) THEN 97 zfjulday = fjulday + rdt / rday 98 IF( ABS(zfjulday - REAL(NINT(zfjulday),wp)) < 0.1 / rday ) zfjulday = REAL(NINT(zfjulday),wp) ! avoid truncation error 99 CALL ju2ymds( zfjulday, iyear, imonth, iday, zsec ) 100 WRITE(clkt, '(i4.4,2i2.2)') iyear, imonth, iday 101 ELSE 102 IF( nitrst > 999999999 ) THEN ; WRITE(clkt, * ) nitrst 103 ELSE ; WRITE(clkt, '(i8.8)') nitrst 104 ENDIF 94 105 ENDIF 95 106 ! create the file … … 173 184 END IF 174 185 ENDIF 186 CALL iom_rstput( kt, nitrst, numrow, 'neos' , REAL(neos) , ldxios = lwxios) ! equation of state 187 !CALL iom_rstput( kt, nitrst, numrow, 'neos' , neos , ktype = jp_i1, ldxios = lwxios) ! equation of state 188 175 189 176 190 IF (ln_diurnal) CALL iom_rstput( kt, nitrst, numrow, 'Dsst', x_dsst, ldxios = lwxios ) … … 249 263 !!---------------------------------------------------------------------- 250 264 REAL(wp) :: zrdt 265 REAL(wp) :: zeos 251 266 INTEGER :: jk 252 267 REAL(wp), DIMENSION(jpi, jpj, jpk) :: w3d … … 255 270 CALL rst_read_open ! open restart for reading (if not already opened) 256 271 272 IF ( ln_rst_eos ) THEN 273 ! Check equation of state used is consistent with the restart 274 IF( iom_varid( numror, 'neos') == -1) THEN 275 CALL ctl_stop( 'restart, rst_read: variable neos not found. STOP check that the equations of state in the restart file and in the namelist nameos are consistent and use ln_rst_eos=F') 276 ELSE 277 CALL iom_get( numror, 'neos', zeos, ldxios = lrxios ) 278 IF ( INT(zeos) /= neos ) CALL ctl_stop( 'restart, rst_read: equation of state used in restart file differs from namelist nameos') 279 ENDIF 280 ENDIF 281 257 282 ! Check dynamics and tracer time-step consistency and force Euler restart if changed 258 IF( iom_varid( numror, 'rdt', ldstop = .FALSE. ) > 0 ) THEN 283 IF( iom_varid( numror, 'rdt', ldstop = .FALSE. ) > 0 ) THEN 259 284 CALL iom_get( numror, 'rdt', zrdt, ldxios = lrxios ) 260 285 IF( zrdt /= rdt ) neuler = 0 -
NEMO/branches/UKMO/NEMO_4.0.4_GO8_package/src/OCE/SBC/sbcssm.F90
r14075 r14078 57 57 REAL(wp) :: zcoef, zf_sbc ! local scalar 58 58 REAL(wp), DIMENSION(jpi,jpj,jpts) :: zts 59 CHARACTER(len=4),SAVE :: stype 59 60 !!--------------------------------------------------------------------- 61 IF( kt == nit000 ) THEN 62 IF( ln_TEOS10 ) THEN 63 stype='abs' ! teos-10: using absolute salinity (sst is converted to potential temperature for the surface module) 64 ELSE IF( ln_EOS80 ) THEN 65 stype='pra' ! eos-80: using practical salinity 66 ELSE IF ( ln_SEOS) THEN 67 stype='seos' ! seos using Simplified Equation of state (sst is converted to potential temperature for the surface module) 68 ENDIF 69 ENDIF 60 70 ! 61 71 ! !* surface T-, U-, V- ocean level variables (T, S, depth, velocity) … … 174 184 CALL iom_put( 'ssu_m', ssu_m ) 175 185 CALL iom_put( 'ssv_m', ssv_m ) 176 CALL iom_put( 'sst_m ', sst_m )177 CALL iom_put( 'sss_m ', sss_m )186 CALL iom_put( 'sst_m_pot', sst_m ) 187 CALL iom_put( 'sss_m_'//stype, sss_m ) 178 188 CALL iom_put( 'ssh_m', ssh_m ) 179 189 CALL iom_put( 'e3t_m', e3t_m ) -
NEMO/branches/UKMO/NEMO_4.0.4_GO8_package/src/OCE/ZDF/zdfmxl.F90
r14075 r14078 15 15 USE trc_oce , ONLY: l_offline ! ocean space and time domain variables 16 16 USE zdf_oce ! ocean vertical physics 17 USE eosbn2 ! for zdf_mxl_zint 17 18 ! 18 19 USE in_out_manager ! I/O manager … … 31 32 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hmlp !: mixed layer depth (rho=rho0+zdcrit) [m] (used by LDF) 32 33 REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: hmlpt !: depth of the last T-point inside the mixed layer [m] (used by LDF) 34 REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: hmld_zint !: vertically-interpolated mixed layer depth [m] 35 REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: htc_mld ! Heat content of hmld_zint 36 LOGICAL, PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: ll_found ! Is T_b to be found by interpolation ? 37 LOGICAL, PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: ll_belowml ! Flag points below mixed layer when ll_found=F 33 38 34 39 REAL(wp), PUBLIC :: rho_c = 0.01_wp !: density criterion for mixed layer depth 35 40 REAL(wp), PUBLIC :: avt_c = 5.e-4_wp ! Kz criterion for the turbocline depth 41 42 TYPE, PUBLIC :: MXL_ZINT !: Structure for MLD defs 43 INTEGER :: mld_type ! mixed layer type 44 REAL(wp) :: zref ! depth of initial T_ref 45 REAL(wp) :: dT_crit ! Critical temp diff 46 REAL(wp) :: iso_frac ! Fraction of rn_dT_crit 47 END TYPE MXL_ZINT 36 48 37 49 !!---------------------------------------------------------------------- … … 48 60 zdf_mxl_alloc = 0 ! set to zero if no array to be allocated 49 61 IF( .NOT. ALLOCATED( nmln ) ) THEN 50 ALLOCATE( nmln(jpi,jpj), hmld(jpi,jpj), hmlp(jpi,jpj), hmlpt(jpi,jpj), STAT= zdf_mxl_alloc ) 62 ALLOCATE( nmln(jpi,jpj), hmld(jpi,jpj), hmlp(jpi,jpj), hmlpt(jpi,jpj), hmld_zint(jpi,jpj), & 63 & htc_mld(jpi,jpj), ll_found(jpi,jpj), ll_belowml(jpi,jpj,jpk), STAT= zdf_mxl_alloc ) 51 64 ! 52 65 CALL mpp_sum ( 'zdfmxl', zdf_mxl_alloc ) … … 137 150 ENDIF 138 151 ! 152 ! Vertically-interpolated mixed-layer depth diagnostic 153 CALL zdf_mxl_zint( kt ) 154 ! 139 155 IF(ln_ctl) CALL prt_ctl( tab2d_1=REAL(nmln,wp), clinfo1=' nmln : ', tab2d_2=hmlp, clinfo2=' hmlp : ' ) 140 156 ! 141 157 END SUBROUTINE zdf_mxl 158 159 SUBROUTINE zdf_mxl_zint_mld( sf ) 160 !!---------------------------------------------------------------------------------- 161 !! *** ROUTINE zdf_mxl_zint_mld *** 162 ! 163 ! Calculate vertically-interpolated mixed layer depth diagnostic. 164 ! 165 ! This routine can calculate the mixed layer depth diagnostic suggested by 166 ! Kara et al, 2000, JGR, 105, 16803, but is more general and can calculate 167 ! vertically-interpolated mixed-layer depth diagnostics with other parameter 168 ! settings set in the namzdf_mldzint namelist. 169 ! 170 ! If mld_type=1 the mixed layer depth is calculated as the depth at which the 171 ! density has increased by an amount equivalent to a temperature difference of 172 ! 0.8C at the surface. 173 ! 174 ! For other values of mld_type the mixed layer is calculated as the depth at 175 ! which the temperature differs by 0.8C from the surface temperature. 176 ! 177 ! David Acreman, Daley Calvert 178 ! 179 !!----------------------------------------------------------------------------------- 180 181 TYPE(MXL_ZINT), INTENT(in) :: sf 182 183 ! Diagnostic criteria 184 INTEGER :: nn_mld_type ! mixed layer type 185 REAL(wp) :: rn_zref ! depth of initial T_ref 186 REAL(wp) :: rn_dT_crit ! Critical temp diff 187 REAL(wp) :: rn_iso_frac ! Fraction of rn_dT_crit used 188 189 ! Local variables 190 REAL(wp), PARAMETER :: zepsilon = 1.e-30 ! local small value 191 INTEGER, DIMENSION(jpi,jpj) :: ikmt ! number of active tracer levels 192 INTEGER, DIMENSION(jpi,jpj) :: ik_ref ! index of reference level 193 INTEGER, DIMENSION(jpi,jpj) :: ik_iso ! index of last uniform temp level 194 REAL, DIMENSION(jpi,jpj,jpk) :: zT ! Temperature or density 195 REAL, DIMENSION(jpi,jpj) :: ppzdep ! depth for use in calculating d(rho) 196 REAL, DIMENSION(jpi,jpj) :: zT_ref ! reference temperature 197 REAL :: zT_b ! base temperature 198 REAL, DIMENSION(jpi,jpj,jpk) :: zdTdz ! gradient of zT 199 REAL, DIMENSION(jpi,jpj,jpk) :: zmoddT ! Absolute temperature difference 200 REAL :: zdz ! depth difference 201 REAL :: zdT ! temperature difference 202 REAL, DIMENSION(jpi,jpj) :: zdelta_T ! difference critereon 203 REAL, DIMENSION(jpi,jpj) :: zRHO1, zRHO2 ! Densities 204 INTEGER :: ji, jj, jk ! loop counter 205 206 !!------------------------------------------------------------------------------------- 207 ! 208 ! Unpack structure 209 nn_mld_type = sf%mld_type 210 rn_zref = sf%zref 211 rn_dT_crit = sf%dT_crit 212 rn_iso_frac = sf%iso_frac 213 214 ! Set the mixed layer depth criterion at each grid point 215 IF( nn_mld_type == 0 ) THEN 216 zdelta_T(:,:) = rn_dT_crit 217 zT(:,:,:) = rhop(:,:,:) 218 ELSE IF( nn_mld_type == 1 ) THEN 219 ppzdep(:,:)=0.0 220 call eos ( tsn(:,:,1,:), ppzdep(:,:), zRHO1(:,:) ) 221 ! Use zT temporarily as a copy of tsn with rn_dT_crit added to SST 222 ! [assumes number of tracers less than number of vertical levels] 223 zT(:,:,1:jpts)=tsn(:,:,1,1:jpts) 224 zT(:,:,jp_tem)=zT(:,:,1)+rn_dT_crit 225 CALL eos( zT(:,:,1:jpts), ppzdep(:,:), zRHO2(:,:) ) 226 zdelta_T(:,:) = abs( zRHO1(:,:) - zRHO2(:,:) ) * rau0 227 ! RHO from eos (2d version) doesn't calculate north or east halo: 228 CALL lbc_lnk( 'zdfmxl', zdelta_T, 'T', 1. ) 229 zT(:,:,:) = rhop(:,:,:) 230 ELSE 231 zdelta_T(:,:) = rn_dT_crit 232 zT(:,:,:) = tsn(:,:,:,jp_tem) 233 END IF 234 235 ! Calculate the gradient of zT and absolute difference for use later 236 DO jk = 1 ,jpk-2 237 zdTdz(:,:,jk) = ( zT(:,:,jk+1) - zT(:,:,jk) ) / e3w_n(:,:,jk+1) 238 zmoddT(:,:,jk) = abs( zT(:,:,jk+1) - zT(:,:,jk) ) 239 END DO 240 241 ! Find density/temperature at the reference level (Kara et al use 10m). 242 ! ik_ref is the index of the box centre immediately above or at the reference level 243 ! Find rn_zref in the array of model level depths and find the ref 244 ! density/temperature by linear interpolation. 245 DO jk = jpkm1, 2, -1 246 WHERE ( gdept_n(:,:,jk) > rn_zref ) 247 ik_ref(:,:) = jk - 1 248 zT_ref(:,:) = zT(:,:,jk-1) + zdTdz(:,:,jk-1) * ( rn_zref - gdept_n(:,:,jk-1) ) 249 END WHERE 250 END DO 251 252 ! If the first grid box centre is below the reference level then use the 253 ! top model level to get zT_ref 254 WHERE ( gdept_n(:,:,1) > rn_zref ) 255 zT_ref = zT(:,:,1) 256 ik_ref = 1 257 END WHERE 258 259 ! The number of active tracer levels is 1 less than the number of active w levels 260 ikmt(:,:) = mbkt(:,:) - 1 261 262 ! Initialize / reset 263 ll_found(:,:) = .false. 264 265 IF ( rn_iso_frac - zepsilon > 0. ) THEN 266 ! Search for a uniform density/temperature region where adjacent levels 267 ! differ by less than rn_iso_frac * deltaT. 268 ! ik_iso is the index of the last level in the uniform layer 269 ! ll_found indicates whether the mixed layer depth can be found by interpolation 270 ik_iso(:,:) = ik_ref(:,:) 271 DO jj = 1, nlcj 272 DO ji = 1, nlci 273 !CDIR NOVECTOR 274 DO jk = ik_ref(ji,jj), ikmt(ji,jj)-1 275 IF ( zmoddT(ji,jj,jk) > ( rn_iso_frac * zdelta_T(ji,jj) ) ) THEN 276 ik_iso(ji,jj) = jk 277 ll_found(ji,jj) = ( zmoddT(ji,jj,jk) > zdelta_T(ji,jj) ) 278 EXIT 279 END IF 280 END DO 281 END DO 282 END DO 283 284 ! Use linear interpolation to find depth of mixed layer base where possible 285 hmld_zint(:,:) = rn_zref 286 DO jj = 1, jpj 287 DO ji = 1, jpi 288 IF (ll_found(ji,jj) .and. tmask(ji,jj,1) == 1.0) THEN 289 zdz = abs( zdelta_T(ji,jj) / zdTdz(ji,jj,ik_iso(ji,jj)) ) 290 hmld_zint(ji,jj) = gdept_n(ji,jj,ik_iso(ji,jj)) + zdz 291 END IF 292 END DO 293 END DO 294 END IF 295 296 ! If ll_found = .false. then calculate MLD using difference of zdelta_T 297 ! from the reference density/temperature 298 299 ! Prevent this section from working on land points 300 WHERE ( tmask(:,:,1) /= 1.0 ) 301 ll_found = .true. 302 END WHERE 303 304 DO jk=1, jpk 305 ll_belowml(:,:,jk) = abs( zT(:,:,jk) - zT_ref(:,:) ) >= zdelta_T(:,:) 306 END DO 307 308 ! Set default value where interpolation cannot be used (ll_found=false) 309 DO jj = 1, jpj 310 DO ji = 1, jpi 311 IF ( .not. ll_found(ji,jj) ) hmld_zint(ji,jj) = gdept_n(ji,jj,ikmt(ji,jj)) 312 END DO 313 END DO 314 315 DO jj = 1, jpj 316 DO ji = 1, jpi 317 !CDIR NOVECTOR 318 DO jk = ik_ref(ji,jj)+1, ikmt(ji,jj) 319 IF ( ll_found(ji,jj) ) EXIT 320 IF ( ll_belowml(ji,jj,jk) ) THEN 321 zT_b = zT_ref(ji,jj) + zdelta_T(ji,jj) * SIGN(1.0, zdTdz(ji,jj,jk-1) ) 322 zdT = zT_b - zT(ji,jj,jk-1) 323 zdz = zdT / zdTdz(ji,jj,jk-1) 324 hmld_zint(ji,jj) = gdept_n(ji,jj,jk-1) + zdz 325 EXIT 326 END IF 327 END DO 328 END DO 329 END DO 330 331 hmld_zint(:,:) = hmld_zint(:,:)*tmask(:,:,1) 332 ! 333 END SUBROUTINE zdf_mxl_zint_mld 334 335 SUBROUTINE zdf_mxl_zint_htc( kt ) 336 !!---------------------------------------------------------------------- 337 !! *** ROUTINE zdf_mxl_zint_htc *** 338 !! 339 !! ** Purpose : 340 !! 341 !! ** Method : 342 !!---------------------------------------------------------------------- 343 344 INTEGER, INTENT(in) :: kt ! ocean time-step index 345 346 INTEGER :: ji, jj, jk 347 INTEGER :: ikmax 348 REAL(wp) :: zc, zcoef 349 ! 350 INTEGER, ALLOCATABLE, DIMENSION(:,:) :: ilevel 351 REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zthick_0, zthick 352 353 !!---------------------------------------------------------------------- 354 355 IF( .NOT. ALLOCATED(ilevel) ) THEN 356 ALLOCATE( ilevel(jpi,jpj), zthick_0(jpi,jpj), & 357 & zthick(jpi,jpj), STAT=ji ) 358 IF( lk_mpp ) CALL mpp_sum( 'zdfmxl', ji ) 359 IF( ji /= 0 ) CALL ctl_stop( 'STOP', 'zdf_mxl_zint_htc : unable to allocate arrays' ) 360 ENDIF 361 362 ! Find last whole model T level above the MLD 363 ilevel(:,:) = 0 364 zthick_0(:,:) = 0._wp 365 366 DO jk = 1, jpkm1 367 DO jj = 1, jpj 368 DO ji = 1, jpi 369 zthick_0(ji,jj) = zthick_0(ji,jj) + e3t_n(ji,jj,jk) 370 IF( zthick_0(ji,jj) < hmld_zint(ji,jj) ) ilevel(ji,jj) = jk 371 END DO 372 END DO 373 WRITE(numout,*) 'zthick_0(jk =',jk,') =',zthick_0(2,2) 374 WRITE(numout,*) 'gdepw_n(jk+1 =',jk+1,') =',gdepw_n(2,2,jk+1) 375 END DO 376 377 ! Surface boundary condition 378 IF( ln_linssh ) THEN ; zthick(:,:) = sshn(:,:) ; htc_mld(:,:) = tsn(:,:,1,jp_tem) * sshn(:,:) * tmask(:,:,1) 379 ELSE ; zthick(:,:) = 0._wp ; htc_mld(:,:) = 0._wp 380 ENDIF 381 382 ! Deepest whole T level above the MLD 383 ikmax = MIN( MAXVAL( ilevel(:,:) ), jpkm1 ) 384 385 ! Integration down to last whole model T level 386 DO jk = 1, ikmax 387 DO jj = 1, jpj 388 DO ji = 1, jpi 389 zc = e3t_n(ji,jj,jk) * REAL( MIN( MAX( 0, ilevel(ji,jj) - jk + 1 ) , 1 ) ) ! 0 below ilevel 390 zthick(ji,jj) = zthick(ji,jj) + zc 391 htc_mld(ji,jj) = htc_mld(ji,jj) + zc * tsn(ji,jj,jk,jp_tem) * tmask(ji,jj,jk) 392 END DO 393 END DO 394 END DO 395 396 ! Subsequent partial T level 397 zthick(:,:) = hmld_zint(:,:) - zthick(:,:) ! remaining thickness to reach MLD 398 399 DO jj = 1, jpj 400 DO ji = 1, jpi 401 htc_mld(ji,jj) = htc_mld(ji,jj) + tsn(ji,jj,ilevel(ji,jj)+1,jp_tem) & 402 & * MIN( e3t_n(ji,jj,ilevel(ji,jj)+1), zthick(ji,jj) ) * tmask(ji,jj,ilevel(ji,jj)+1) 403 END DO 404 END DO 405 406 WRITE(numout,*) 'htc_mld(after) =',htc_mld(2,2) 407 408 ! Convert to heat content 409 zcoef = rau0 * rcp 410 htc_mld(:,:) = zcoef * htc_mld(:,:) 411 412 END SUBROUTINE zdf_mxl_zint_htc 413 414 SUBROUTINE zdf_mxl_zint( kt ) 415 !!---------------------------------------------------------------------- 416 !! *** ROUTINE zdf_mxl_zint *** 417 !! 418 !! ** Purpose : 419 !! 420 !! ** Method : 421 !!---------------------------------------------------------------------- 422 423 INTEGER, INTENT(in) :: kt ! ocean time-step index 424 425 INTEGER :: ios 426 INTEGER :: jn 427 428 INTEGER :: nn_mld_diag = 0 ! number of diagnostics 429 430 CHARACTER(len=1) :: cmld 431 432 TYPE(MXL_ZINT) :: sn_mld1, sn_mld2, sn_mld3, sn_mld4, sn_mld5 433 TYPE(MXL_ZINT), SAVE, DIMENSION(5) :: mld_diags 434 435 NAMELIST/namzdf_mldzint/ nn_mld_diag, sn_mld1, sn_mld2, sn_mld3, sn_mld4, sn_mld5 436 437 !!---------------------------------------------------------------------- 438 439 IF( kt == nit000 ) THEN 440 REWIND( numnam_ref ) ! Namelist namzdf_mldzint in reference namelist 441 READ ( numnam_ref, namzdf_mldzint, IOSTAT = ios, ERR = 901) 442 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_mldzint in reference namelist' ) 443 444 REWIND( numnam_cfg ) ! Namelist namzdf_mldzint in configuration namelist 445 READ ( numnam_cfg, namzdf_mldzint, IOSTAT = ios, ERR = 902 ) 446 902 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namzdf_mldzint in configuration namelist' ) 447 IF(lwm) WRITE ( numond, namzdf_mldzint ) 448 449 IF( nn_mld_diag > 5 ) CALL ctl_stop( 'STOP', 'zdf_mxl_ini: Specify no more than 5 MLD definitions' ) 450 451 mld_diags(1) = sn_mld1 452 mld_diags(2) = sn_mld2 453 mld_diags(3) = sn_mld3 454 mld_diags(4) = sn_mld4 455 mld_diags(5) = sn_mld5 456 457 IF( lwp .AND. (nn_mld_diag > 0) ) THEN 458 WRITE(numout,*) '=============== Vertically-interpolated mixed layer ================' 459 WRITE(numout,*) '(Diagnostic number, nn_mld_type, rn_zref, rn_dT_crit, rn_iso_frac)' 460 DO jn = 1, nn_mld_diag 461 WRITE(numout,*) 'MLD criterion',jn,':' 462 WRITE(numout,*) ' nn_mld_type =', mld_diags(jn)%mld_type 463 WRITE(numout,*) ' rn_zref =' , mld_diags(jn)%zref 464 WRITE(numout,*) ' rn_dT_crit =' , mld_diags(jn)%dT_crit 465 WRITE(numout,*) ' rn_iso_frac =', mld_diags(jn)%iso_frac 466 END DO 467 WRITE(numout,*) '====================================================================' 468 ENDIF 469 ENDIF 470 471 IF( nn_mld_diag > 0 ) THEN 472 DO jn = 1, nn_mld_diag 473 WRITE(cmld,'(I1)') jn 474 IF( iom_use( "mldzint_"//cmld ) .OR. iom_use( "mldhtc_"//cmld ) ) THEN 475 CALL zdf_mxl_zint_mld( mld_diags(jn) ) 476 477 IF( iom_use( "mldzint_"//cmld ) ) THEN 478 CALL iom_put( "mldzint_"//cmld, hmld_zint(:,:) ) 479 ENDIF 480 481 IF( iom_use( "mldhtc_"//cmld ) ) THEN 482 CALL zdf_mxl_zint_htc( kt ) 483 CALL iom_put( "mldhtc_"//cmld , htc_mld(:,:) ) 484 ENDIF 485 ENDIF 486 END DO 487 ENDIF 488 489 END SUBROUTINE zdf_mxl_zint 142 490 143 491 !!====================================================================== -
NEMO/branches/UKMO/NEMO_4.0.4_GO8_package/src/OCE/step.F90
r14075 r14078 207 207 IF( ln_diaptr ) CALL dia_ptr ! Poleward adv/ldf TRansports diagnostics 208 208 IF( ln_diaharm ) CALL dia_harm( kstp ) ! Tidal harmonic analysis 209 CALL dia_prod( kstp ) ! ocean model: product diagnostics 209 210 CALL dia_wri ( kstp ) ! ocean model: outputs 210 211 ! -
NEMO/branches/UKMO/NEMO_4.0.4_GO8_package/src/OCE/step_oce.F90
r14075 r14078 80 80 USE diahsb ! heat, salt and volume budgets (dia_hsb routine) 81 81 USE diaharm 82 USE diaprod 82 83 USE diacfl 83 84 USE diaobs ! Observation operator -
NEMO/branches/UKMO/NEMO_4.0.4_GO8_package/src/SAS/sbcssm.F90
r14075 r14078 77 77 REAL(wp) :: ztinta ! ratio applied to after records when doing time interpolation 78 78 REAL(wp) :: ztintb ! ratio applied to before records when doing time interpolation 79 !!---------------------------------------------------------------------- 79 CHARACTER(len=4),SAVE :: stype 80 !!--------------------------------------------------------------------- 81 IF( kt == nit000 ) THEN 82 IF( ln_TEOS10 ) THEN 83 stype='abs' ! teos-10: using absolute salinity (sst is converted to potential temperature for the surface module) 84 ELSE IF( ln_EOS80 ) THEN 85 stype='pra' ! eos-80: using practical salinity 86 ELSE IF ( ln_SEOS) THEN 87 stype='seos' ! seos using Simplified Equation of state (sst is converted to potential temperature for the surface module) 88 ENDIF 89 ENDIF 80 90 ! 81 91 IF( ln_timing ) CALL timing_start( 'sbc_ssm') … … 144 154 CALL iom_put( 'ssu_m', ssu_m ) 145 155 CALL iom_put( 'ssv_m', ssv_m ) 146 CALL iom_put( 'sst_m ', sst_m )147 CALL iom_put( 'sss_m ', sss_m )156 CALL iom_put( 'sst_m_pot', sst_m ) 157 CALL iom_put( 'sss_m_'//stype, sss_m ) 148 158 CALL iom_put( 'ssh_m', ssh_m ) 149 159 IF( .NOT.ln_linssh ) CALL iom_put( 'e3t_m', e3t_m )
Note: See TracChangeset
for help on using the changeset viewer.