Changeset 5965 for branches/2014/dev_r4650_UKMO14.5_SST_BIAS_CORRECTION/NEMOGCM/NEMO/OPA_SRC/SBC/sbcice_lim.F90
- Timestamp:
- 2015-12-01T16:35:30+01:00 (8 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/2014/dev_r4650_UKMO14.5_SST_BIAS_CORRECTION/NEMOGCM/NEMO/OPA_SRC/SBC/sbcice_lim.F90
r4333 r5965 12 12 !! 3.4 ! 2011-01 (A Porter) dynamical allocation 13 13 !! - ! 2012-10 (C. Rousset) add lim_diahsb 14 !! 3.6 ! 2014-07 (M. Vancoppenolle, G. Madec, O. Marti) revise coupled interface 14 15 !!---------------------------------------------------------------------- 15 16 #if defined key_lim3 … … 18 19 !!---------------------------------------------------------------------- 19 20 !! sbc_ice_lim : sea-ice model time-stepping and update ocean sbc over ice-covered area 20 !! lim_ctl : alerts in case of ice model crash21 !! lim_prt_state : ice control print at a given grid point22 21 !!---------------------------------------------------------------------- 23 22 USE oce ! ocean dynamics and tracers 24 23 USE dom_oce ! ocean space and time domain 25 USE par_ice ! sea-ice parameters26 24 USE ice ! LIM-3: ice variables 27 USE iceini ! LIM-3: ice initialisation25 USE thd_ice ! LIM-3: thermodynamical variables 28 26 USE dom_ice ! LIM-3: ice domain 29 27 … … 39 37 USE limdyn ! Ice dynamics 40 38 USE limtrp ! Ice transport 39 USE limhdf ! Ice horizontal diffusion 41 40 USE limthd ! Ice thermodynamics 42 USE limitd_th ! Thermodynamics on ice thickness distribution43 41 USE limitd_me ! Mechanics on ice thickness distribution 44 42 USE limsbc ! sea surface boundary condition … … 46 44 USE limwri ! Ice outputs 47 45 USE limrst ! Ice restarts 48 USE limupdate1 49 USE limupdate2 46 USE limupdate1 ! update of global variables 47 USE limupdate2 ! update of global variables 50 48 USE limvar ! Ice variables switch 49 50 USE limmsh ! LIM mesh 51 USE limistate ! LIM initial state 52 USE limthd_sal ! LIM ice thermodynamics: salinity 51 53 52 54 USE c1d ! 1D vertical configuration … … 59 61 USE prtctl ! Print control 60 62 USE lib_fortran ! 63 USE limctl 61 64 62 65 #if defined key_bdy … … 68 71 69 72 PUBLIC sbc_ice_lim ! routine called by sbcmod.F90 73 PUBLIC sbc_lim_init ! routine called by sbcmod.F90 70 74 71 75 !! * Substitutions … … 79 83 CONTAINS 80 84 81 FUNCTION fice_cell_ave ( ptab) 85 !!====================================================================== 86 87 SUBROUTINE sbc_ice_lim( kt, kblk ) 88 !!--------------------------------------------------------------------- 89 !! *** ROUTINE sbc_ice_lim *** 90 !! 91 !! ** Purpose : update the ocean surface boundary condition via the 92 !! Louvain la Neuve Sea Ice Model time stepping 93 !! 94 !! ** Method : ice model time stepping 95 !! - call the ice dynamics routine 96 !! - call the ice advection/diffusion routine 97 !! - call the ice thermodynamics routine 98 !! - call the routine that computes mass and 99 !! heat fluxes at the ice/ocean interface 100 !! - save the outputs 101 !! - save the outputs for restart when necessary 102 !! 103 !! ** Action : - time evolution of the LIM sea-ice model 104 !! - update all sbc variables below sea-ice: 105 !! utau, vtau, taum, wndm, qns , qsr, emp , sfx 106 !!--------------------------------------------------------------------- 107 INTEGER, INTENT(in) :: kt ! ocean time step 108 INTEGER, INTENT(in) :: kblk ! type of bulk (=3 CLIO, =4 CORE, =5 COUPLED) 109 !! 110 INTEGER :: jl ! dummy loop index 111 REAL(wp), POINTER, DIMENSION(:,:,:) :: zalb_os, zalb_cs ! ice albedo under overcast/clear sky 112 REAL(wp), POINTER, DIMENSION(:,:,:) :: zalb_ice ! mean ice albedo (for coupled) 113 REAL(wp), POINTER, DIMENSION(:,: ) :: zutau_ice, zvtau_ice 114 !!---------------------------------------------------------------------- 115 116 IF( nn_timing == 1 ) CALL timing_start('sbc_ice_lim') 117 118 !-----------------------! 119 ! --- Ice time step --- ! 120 !-----------------------! 121 IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN 122 123 ! mean surface ocean current at ice velocity point (C-grid dynamics : U- & V-points as the ocean) 124 u_oce(:,:) = ssu_m(:,:) * umask(:,:,1) 125 v_oce(:,:) = ssv_m(:,:) * vmask(:,:,1) 126 127 ! masked sea surface freezing temperature [Kelvin] (set to rt0 over land) 128 t_bo(:,:) = ( eos_fzp( sss_m ) + rt0 ) * tmask(:,:,1) + rt0 * ( 1._wp - tmask(:,:,1) ) 129 130 ! Mask sea ice surface temperature (set to rt0 over land) 131 DO jl = 1, jpl 132 t_su(:,:,jl) = t_su(:,:,jl) * tmask(:,:,1) + rt0 * ( 1._wp - tmask(:,:,1) ) 133 END DO 134 ! 135 !------------------------------------------------! 136 ! --- Dynamical coupling with the atmosphere --- ! 137 !------------------------------------------------! 138 ! It provides the following fields: 139 ! utau_ice, vtau_ice : surface ice stress (U- & V-points) [N/m2] 140 !----------------------------------------------------------------- 141 SELECT CASE( kblk ) 142 CASE( jp_clio ) ; CALL blk_ice_clio_tau ! CLIO bulk formulation 143 CASE( jp_core ) ; CALL blk_ice_core_tau ! CORE bulk formulation 144 CASE( jp_purecpl ) ; CALL sbc_cpl_ice_tau( utau_ice , vtau_ice ) ! Coupled formulation 145 END SELECT 146 147 IF( ln_mixcpl) THEN ! Case of a mixed Bulk/Coupled formulation 148 CALL wrk_alloc( jpi,jpj , zutau_ice, zvtau_ice) 149 CALL sbc_cpl_ice_tau( zutau_ice , zvtau_ice ) 150 utau_ice(:,:) = utau_ice(:,:) * xcplmask(:,:,0) + zutau_ice(:,:) * ( 1. - xcplmask(:,:,0) ) 151 vtau_ice(:,:) = vtau_ice(:,:) * xcplmask(:,:,0) + zvtau_ice(:,:) * ( 1. - xcplmask(:,:,0) ) 152 CALL wrk_dealloc( jpi,jpj , zutau_ice, zvtau_ice) 153 ENDIF 154 155 !-------------------------------------------------------! 156 ! --- ice dynamics and transport (except in 1D case) ---! 157 !-------------------------------------------------------! 158 numit = numit + nn_fsbc ! Ice model time step 159 ! 160 CALL sbc_lim_bef ! Store previous ice values 161 CALL sbc_lim_diag0 ! set diag of mass, heat and salt fluxes to 0 162 CALL lim_rst_opn( kt ) ! Open Ice restart file 163 ! 164 IF( .NOT. lk_c1d ) THEN 165 ! 166 CALL lim_dyn( kt ) ! Ice dynamics ( rheology/dynamics ) 167 ! 168 CALL lim_trp( kt ) ! Ice transport ( Advection/diffusion ) 169 ! 170 IF( nn_monocat /= 2 ) CALL lim_itd_me ! Mechanical redistribution ! (ridging/rafting) 171 ! 172 #if defined key_bdy 173 CALL bdy_ice_lim( kt ) ! bdy ice thermo 174 IF( ln_icectl ) CALL lim_prt( kt, iiceprt, jiceprt, 1, ' - ice thermo bdy - ' ) 175 #endif 176 ! 177 CALL lim_update1( kt ) ! Corrections 178 ! 179 ENDIF 180 181 ! previous lead fraction and ice volume for flux calculations 182 CALL sbc_lim_bef 183 CALL lim_var_glo2eqv ! ht_i and ht_s for ice albedo calculation 184 CALL lim_var_agg(1) ! at_i for coupling (via pfrld) 185 pfrld(:,:) = 1._wp - at_i(:,:) 186 phicif(:,:) = vt_i(:,:) 187 188 !------------------------------------------------------! 189 ! --- Thermodynamical coupling with the atmosphere --- ! 190 !------------------------------------------------------! 191 ! It provides the following fields: 192 ! qsr_ice , qns_ice : solar & non solar heat flux over ice (T-point) [W/m2] 193 ! qla_ice : latent heat flux over ice (T-point) [W/m2] 194 ! dqns_ice, dqla_ice : non solar & latent heat sensistivity (T-point) [W/m2] 195 ! tprecip , sprecip : total & solid precipitation (T-point) [Kg/m2/s] 196 ! fr1_i0 , fr2_i0 : 1sr & 2nd fraction of qsr penetration in ice [%] 197 !---------------------------------------------------------------------------------------- 198 CALL wrk_alloc( jpi,jpj,jpl, zalb_os, zalb_cs, zalb_ice ) 199 CALL albedo_ice( t_su, ht_i, ht_s, zalb_cs, zalb_os ) ! cloud-sky and overcast-sky ice albedos 200 201 SELECT CASE( kblk ) 202 CASE( jp_clio ) ! CLIO bulk formulation 203 ! In CLIO the cloud fraction is read in the climatology and the all-sky albedo 204 ! (zalb_ice) is computed within the bulk routine 205 CALL blk_ice_clio_flx( t_su, zalb_cs, zalb_os, zalb_ice ) 206 IF( ln_mixcpl ) CALL sbc_cpl_ice_flx( p_frld=pfrld, palbi=zalb_ice, psst=sst_m, pist=t_su ) 207 IF( nn_limflx /= 2 ) CALL ice_lim_flx( t_su, zalb_ice, qns_ice, qsr_ice, dqns_ice, evap_ice, devap_ice, nn_limflx ) 208 CASE( jp_core ) ! CORE bulk formulation 209 ! albedo depends on cloud fraction because of non-linear spectral effects 210 zalb_ice(:,:,:) = ( 1. - cldf_ice ) * zalb_cs(:,:,:) + cldf_ice * zalb_os(:,:,:) 211 CALL blk_ice_core_flx( t_su, zalb_ice ) 212 IF( ln_mixcpl ) CALL sbc_cpl_ice_flx( p_frld=pfrld, palbi=zalb_ice, psst=sst_m, pist=t_su ) 213 IF( nn_limflx /= 2 ) CALL ice_lim_flx( t_su, zalb_ice, qns_ice, qsr_ice, dqns_ice, evap_ice, devap_ice, nn_limflx ) 214 CASE ( jp_purecpl ) 215 ! albedo depends on cloud fraction because of non-linear spectral effects 216 zalb_ice(:,:,:) = ( 1. - cldf_ice ) * zalb_cs(:,:,:) + cldf_ice * zalb_os(:,:,:) 217 CALL sbc_cpl_ice_flx( p_frld=pfrld, palbi=zalb_ice, psst=sst_m, pist=t_su ) 218 ! clem: evap_ice is forced to 0 in coupled mode for now 219 ! but it needs to be changed (along with modif in limthd_dh) once heat flux from evap will be avail. from atm. models 220 evap_ice (:,:,:) = 0._wp ; devap_ice (:,:,:) = 0._wp 221 IF( nn_limflx == 2 ) CALL ice_lim_flx( t_su, zalb_ice, qns_ice, qsr_ice, dqns_ice, evap_ice, devap_ice, nn_limflx ) 222 END SELECT 223 CALL wrk_dealloc( jpi,jpj,jpl, zalb_os, zalb_cs, zalb_ice ) 224 225 !----------------------------! 226 ! --- ice thermodynamics --- ! 227 !----------------------------! 228 CALL lim_thd( kt ) ! Ice thermodynamics 229 ! 230 CALL lim_update2( kt ) ! Corrections 231 ! 232 CALL lim_sbc_flx( kt ) ! Update surface ocean mass, heat and salt fluxes 233 ! 234 IF(ln_limdiaout) CALL lim_diahsb ! Diagnostics and outputs 235 ! 236 CALL lim_wri( 1 ) ! Ice outputs 237 ! 238 IF( kt == nit000 .AND. ln_rstart ) & 239 & CALL iom_close( numrir ) ! close input ice restart file 240 ! 241 IF( lrst_ice ) CALL lim_rst_write( kt ) ! Ice restart file 242 ! 243 IF( ln_icectl ) CALL lim_ctl( kt ) ! alerts in case of model crash 244 ! 245 ENDIF ! End sea-ice time step only 246 247 !-------------------------! 248 ! --- Ocean time step --- ! 249 !-------------------------! 250 ! Update surface ocean stresses (only in ice-dynamic case) otherwise the atm.-ocean stresses are used everywhere 251 IF( ln_limdyn ) CALL lim_sbc_tau( kt, ub(:,:,1), vb(:,:,1) ) ! using before instantaneous surf. currents 252 !!gm remark, the ocean-ice stress is not saved in ice diag call above ..... find a solution!!! 253 ! 254 IF( nn_timing == 1 ) CALL timing_stop('sbc_ice_lim') 255 ! 256 END SUBROUTINE sbc_ice_lim 257 258 259 SUBROUTINE sbc_lim_init 260 !!---------------------------------------------------------------------- 261 !! *** ROUTINE sbc_lim_init *** 262 !! 263 !! ** purpose : Allocate all the dynamic arrays of the LIM-3 modules 264 !!---------------------------------------------------------------------- 265 INTEGER :: ierr 266 !!---------------------------------------------------------------------- 267 IF(lwp) WRITE(numout,*) 268 IF(lwp) WRITE(numout,*) 'sbc_ice_lim : update ocean surface boudary condition' 269 IF(lwp) WRITE(numout,*) '~~~~~~~~~~~ via Louvain la Neuve Ice Model (LIM-3) time stepping' 270 ! 271 ! Open the reference and configuration namelist files and namelist output file 272 CALL ctl_opn( numnam_ice_ref, 'namelist_ice_ref', 'OLD', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp ) 273 CALL ctl_opn( numnam_ice_cfg, 'namelist_ice_cfg', 'OLD', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp ) 274 IF(lwm) CALL ctl_opn( numoni, 'output.namelist.ice', 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp, 1 ) 275 276 CALL ice_run ! set some ice run parameters 277 ! 278 ! ! Allocate the ice arrays 279 ierr = ice_alloc () ! ice variables 280 ierr = ierr + dom_ice_alloc () ! domain 281 ierr = ierr + sbc_ice_alloc () ! surface forcing 282 ierr = ierr + thd_ice_alloc () ! thermodynamics 283 ierr = ierr + lim_itd_me_alloc () ! ice thickness distribution - mechanics 284 ! 285 IF( lk_mpp ) CALL mpp_sum( ierr ) 286 IF( ierr /= 0 ) CALL ctl_stop('STOP', 'sbc_lim_init : unable to allocate ice arrays') 287 ! 288 ! ! adequation jpk versus ice/snow layers/categories 289 IF( jpl > jpk .OR. (nlay_i+1) > jpk .OR. nlay_s > jpk ) & 290 & CALL ctl_stop( 'STOP', & 291 & 'sbc_lim_init: the 3rd dimension of workspace arrays is too small.', & 292 & 'use more ocean levels or less ice/snow layers/categories.' ) 293 ! 294 CALL lim_itd_init ! ice thickness distribution initialization 295 ! 296 CALL lim_hdf_init ! set ice horizontal diffusion computation parameters 297 ! 298 CALL lim_thd_init ! set ice thermodynics parameters 299 ! 300 CALL lim_thd_sal_init ! set ice salinity parameters 301 ! 302 CALL lim_msh ! ice mesh initialization 303 ! 304 CALL lim_itd_me_init ! ice thickness distribution initialization for mecanical deformation 305 ! ! Initial sea-ice state 306 IF( .NOT. ln_rstart ) THEN ! start from rest: sea-ice deduced from sst 307 numit = 0 308 numit = nit000 - 1 309 CALL lim_istate 310 ELSE ! start from a restart file 311 CALL lim_rst_read 312 numit = nit000 - 1 313 ENDIF 314 CALL lim_var_agg(1) 315 CALL lim_var_glo2eqv 316 ! 317 CALL lim_sbc_init ! ice surface boundary condition 318 ! 319 fr_i(:,:) = at_i(:,:) ! initialisation of sea-ice fraction 320 tn_ice(:,:,:) = t_su(:,:,:) ! initialisation of surface temp for coupled simu 321 ! 322 nstart = numit + nn_fsbc 323 nitrun = nitend - nit000 + 1 324 nlast = numit + nitrun 325 ! 326 IF( nstock == 0 ) nstock = nlast + 1 327 ! 328 END SUBROUTINE sbc_lim_init 329 330 331 SUBROUTINE ice_run 332 !!------------------------------------------------------------------- 333 !! *** ROUTINE ice_run *** 334 !! 335 !! ** Purpose : Definition some run parameter for ice model 336 !! 337 !! ** Method : Read the namicerun namelist and check the parameter 338 !! values called at the first timestep (nit000) 339 !! 340 !! ** input : Namelist namicerun 341 !!------------------------------------------------------------------- 342 INTEGER :: ios ! Local integer output status for namelist read 343 NAMELIST/namicerun/ jpl, nlay_i, nlay_s, cn_icerst_in, cn_icerst_indir, cn_icerst_out, cn_icerst_outdir, & 344 & ln_limdyn, rn_amax, ln_limdiahsb, ln_limdiaout, ln_icectl, iiceprt, jiceprt 345 !!------------------------------------------------------------------- 346 ! 347 REWIND( numnam_ice_ref ) ! Namelist namicerun in reference namelist : Parameters for ice 348 READ ( numnam_ice_ref, namicerun, IOSTAT = ios, ERR = 901) 349 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicerun in reference namelist', lwp ) 350 351 REWIND( numnam_ice_cfg ) ! Namelist namicerun in configuration namelist : Parameters for ice 352 READ ( numnam_ice_cfg, namicerun, IOSTAT = ios, ERR = 902 ) 353 902 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicerun in configuration namelist', lwp ) 354 IF(lwm) WRITE ( numoni, namicerun ) 355 ! 356 ! 357 IF(lwp) THEN ! control print 358 WRITE(numout,*) 359 WRITE(numout,*) 'ice_run : ice share parameters for dynamics/advection/thermo of sea-ice' 360 WRITE(numout,*) ' ~~~~~~' 361 WRITE(numout,*) ' number of ice categories = ', jpl 362 WRITE(numout,*) ' number of ice layers = ', nlay_i 363 WRITE(numout,*) ' number of snow layers = ', nlay_s 364 WRITE(numout,*) ' switch for ice dynamics (1) or not (0) ln_limdyn = ', ln_limdyn 365 WRITE(numout,*) ' maximum ice concentration = ', rn_amax 366 WRITE(numout,*) ' Diagnose heat/salt budget or not ln_limdiahsb = ', ln_limdiahsb 367 WRITE(numout,*) ' Output heat/salt budget or not ln_limdiaout = ', ln_limdiaout 368 WRITE(numout,*) ' control prints in ocean.out for (i,j)=(iiceprt,jiceprt) = ', ln_icectl 369 WRITE(numout,*) ' i-index for control prints (ln_icectl=true) = ', iiceprt 370 WRITE(numout,*) ' j-index for control prints (ln_icectl=true) = ', jiceprt 371 ENDIF 372 ! 373 ! sea-ice timestep and inverse 374 rdt_ice = nn_fsbc * rdttra(1) 375 r1_rdtice = 1._wp / rdt_ice 376 377 ! inverse of nlay_i and nlay_s 378 r1_nlay_i = 1._wp / REAL( nlay_i, wp ) 379 r1_nlay_s = 1._wp / REAL( nlay_s, wp ) 380 ! 381 #if defined key_bdy 382 IF( lwp .AND. ln_limdiahsb ) CALL ctl_warn('online conservation check activated but it does not work with BDY') 383 #endif 384 ! 385 END SUBROUTINE ice_run 386 387 388 SUBROUTINE lim_itd_init 389 !!------------------------------------------------------------------ 390 !! *** ROUTINE lim_itd_init *** 391 !! 392 !! ** Purpose : Initializes the ice thickness distribution 393 !! ** Method : ... 394 !! ** input : Namelist namiceitd 395 !!------------------------------------------------------------------- 396 INTEGER :: ios ! Local integer output status for namelist read 397 NAMELIST/namiceitd/ nn_catbnd, rn_himean 398 ! 399 INTEGER :: jl ! dummy loop index 400 REAL(wp) :: zc1, zc2, zc3, zx1 ! local scalars 401 REAL(wp) :: zhmax, znum, zden, zalpha ! 402 !!------------------------------------------------------------------ 403 ! 404 REWIND( numnam_ice_ref ) ! Namelist namiceitd in reference namelist : Parameters for ice 405 READ ( numnam_ice_ref, namiceitd, IOSTAT = ios, ERR = 903) 406 903 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namiceitd in reference namelist', lwp ) 407 408 REWIND( numnam_ice_cfg ) ! Namelist namiceitd in configuration namelist : Parameters for ice 409 READ ( numnam_ice_cfg, namiceitd, IOSTAT = ios, ERR = 904 ) 410 904 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namiceitd in configuration namelist', lwp ) 411 IF(lwm) WRITE ( numoni, namiceitd ) 412 ! 413 ! 414 IF(lwp) THEN ! control print 415 WRITE(numout,*) 416 WRITE(numout,*) 'ice_itd : ice cat distribution' 417 WRITE(numout,*) ' ~~~~~~' 418 WRITE(numout,*) ' shape of ice categories distribution nn_catbnd = ', nn_catbnd 419 WRITE(numout,*) ' mean ice thickness in the domain (only active if nn_catbnd=2) rn_himean = ', rn_himean 420 ENDIF 421 422 !---------------------------------- 423 !- Thickness categories boundaries 424 !---------------------------------- 425 IF(lwp) WRITE(numout,*) 426 IF(lwp) WRITE(numout,*) 'lim_itd_init : Initialization of ice cat distribution ' 427 IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~' 428 429 hi_max(:) = 0._wp 430 431 SELECT CASE ( nn_catbnd ) 432 !---------------------- 433 CASE (1) ! tanh function (CICE) 434 !---------------------- 435 zc1 = 3._wp / REAL( jpl, wp ) 436 zc2 = 10._wp * zc1 437 zc3 = 3._wp 438 439 DO jl = 1, jpl 440 zx1 = REAL( jl-1, wp ) / REAL( jpl, wp ) 441 hi_max(jl) = hi_max(jl-1) + zc1 + zc2 * (1._wp + TANH( zc3 * (zx1 - 1._wp ) ) ) 442 END DO 443 444 !---------------------- 445 CASE (2) ! h^(-alpha) function 446 !---------------------- 447 zalpha = 0.05 ! exponent of the transform function 448 449 zhmax = 3.*rn_himean 450 451 DO jl = 1, jpl 452 znum = jpl * ( zhmax+1 )**zalpha 453 zden = ( jpl - jl ) * ( zhmax+1 )**zalpha + jl 454 hi_max(jl) = ( znum / zden )**(1./zalpha) - 1 455 END DO 456 457 END SELECT 458 459 DO jl = 1, jpl 460 hi_mean(jl) = ( hi_max(jl) + hi_max(jl-1) ) * 0.5_wp 461 END DO 462 463 ! Set hi_max(jpl) to a big value to ensure that all ice is thinner than hi_max(jpl) 464 hi_max(jpl) = 99._wp 465 466 IF(lwp) WRITE(numout,*) ' Thickness category boundaries ' 467 IF(lwp) WRITE(numout,*) ' hi_max ', hi_max(0:jpl) 468 ! 469 END SUBROUTINE lim_itd_init 470 471 472 SUBROUTINE ice_lim_flx( ptn_ice, palb_ice, pqns_ice, pqsr_ice, pdqn_ice, pevap_ice, pdevap_ice, k_limflx ) 473 !!--------------------------------------------------------------------- 474 !! *** ROUTINE ice_lim_flx *** 475 !! 476 !! ** Purpose : update the ice surface boundary condition by averaging and / or 477 !! redistributing fluxes on ice categories 478 !! 479 !! ** Method : average then redistribute 480 !! 481 !! ** Action : 482 !!--------------------------------------------------------------------- 483 INTEGER , INTENT(in ) :: k_limflx ! =-1 do nothing; =0 average ; 484 ! =1 average and redistribute ; =2 redistribute 485 REAL(wp), DIMENSION(:,:,:), INTENT(in ) :: ptn_ice ! ice surface temperature 486 REAL(wp), DIMENSION(:,:,:), INTENT(in ) :: palb_ice ! ice albedo 487 REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pqns_ice ! non solar flux 488 REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pqsr_ice ! net solar flux 489 REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pdqn_ice ! non solar flux sensitivity 490 REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pevap_ice ! sublimation 491 REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pdevap_ice ! sublimation sensitivity 492 ! 493 INTEGER :: jl ! dummy loop index 494 ! 495 REAL(wp), POINTER, DIMENSION(:,:) :: zalb_m ! Mean albedo over all categories 496 REAL(wp), POINTER, DIMENSION(:,:) :: ztem_m ! Mean temperature over all categories 497 ! 498 REAL(wp), POINTER, DIMENSION(:,:) :: z_qsr_m ! Mean solar heat flux over all categories 499 REAL(wp), POINTER, DIMENSION(:,:) :: z_qns_m ! Mean non solar heat flux over all categories 500 REAL(wp), POINTER, DIMENSION(:,:) :: z_evap_m ! Mean sublimation over all categories 501 REAL(wp), POINTER, DIMENSION(:,:) :: z_dqn_m ! Mean d(qns)/dT over all categories 502 REAL(wp), POINTER, DIMENSION(:,:) :: z_devap_m ! Mean d(evap)/dT over all categories 503 !!---------------------------------------------------------------------- 504 505 IF( nn_timing == 1 ) CALL timing_start('ice_lim_flx') 506 ! 507 ! 508 SELECT CASE( k_limflx ) !== averaged on all ice categories ==! 509 CASE( 0 , 1 ) 510 CALL wrk_alloc( jpi,jpj, z_qsr_m, z_qns_m, z_evap_m, z_dqn_m, z_devap_m) 511 ! 512 z_qns_m (:,:) = fice_ice_ave ( pqns_ice (:,:,:) ) 513 z_qsr_m (:,:) = fice_ice_ave ( pqsr_ice (:,:,:) ) 514 z_dqn_m (:,:) = fice_ice_ave ( pdqn_ice (:,:,:) ) 515 z_evap_m (:,:) = fice_ice_ave ( pevap_ice (:,:,:) ) 516 z_devap_m(:,:) = fice_ice_ave ( pdevap_ice (:,:,:) ) 517 DO jl = 1, jpl 518 pdqn_ice (:,:,jl) = z_dqn_m(:,:) 519 pdevap_ice(:,:,jl) = z_devap_m(:,:) 520 END DO 521 ! 522 DO jl = 1, jpl 523 pqns_ice (:,:,jl) = z_qns_m(:,:) 524 pqsr_ice (:,:,jl) = z_qsr_m(:,:) 525 pevap_ice(:,:,jl) = z_evap_m(:,:) 526 END DO 527 ! 528 CALL wrk_dealloc( jpi,jpj, z_qsr_m, z_qns_m, z_evap_m, z_dqn_m, z_devap_m) 529 END SELECT 530 531 SELECT CASE( k_limflx ) !== redistribution on all ice categories ==! 532 CASE( 1 , 2 ) 533 CALL wrk_alloc( jpi,jpj, zalb_m, ztem_m ) 534 ! 535 zalb_m(:,:) = fice_ice_ave ( palb_ice (:,:,:) ) 536 ztem_m(:,:) = fice_ice_ave ( ptn_ice (:,:,:) ) 537 DO jl = 1, jpl 538 pqns_ice (:,:,jl) = pqns_ice (:,:,jl) + pdqn_ice (:,:,jl) * ( ptn_ice(:,:,jl) - ztem_m(:,:) ) 539 pevap_ice(:,:,jl) = pevap_ice(:,:,jl) + pdevap_ice(:,:,jl) * ( ptn_ice(:,:,jl) - ztem_m(:,:) ) 540 pqsr_ice (:,:,jl) = pqsr_ice (:,:,jl) * ( 1._wp - palb_ice(:,:,jl) ) / ( 1._wp - zalb_m(:,:) ) 541 END DO 542 ! 543 CALL wrk_dealloc( jpi,jpj, zalb_m, ztem_m ) 544 END SELECT 545 ! 546 IF( nn_timing == 1 ) CALL timing_stop('ice_lim_flx') 547 ! 548 END SUBROUTINE ice_lim_flx 549 550 SUBROUTINE sbc_lim_bef 551 !!---------------------------------------------------------------------- 552 !! *** ROUTINE sbc_lim_bef *** 553 !! 554 !! ** purpose : store ice variables at "before" time step 555 !!---------------------------------------------------------------------- 556 a_i_b (:,:,:) = a_i (:,:,:) ! ice area 557 e_i_b (:,:,:,:) = e_i (:,:,:,:) ! ice thermal energy 558 v_i_b (:,:,:) = v_i (:,:,:) ! ice volume 559 v_s_b (:,:,:) = v_s (:,:,:) ! snow volume 560 e_s_b (:,:,:,:) = e_s (:,:,:,:) ! snow thermal energy 561 smv_i_b(:,:,:) = smv_i(:,:,:) ! salt content 562 oa_i_b (:,:,:) = oa_i (:,:,:) ! areal age content 563 u_ice_b(:,:) = u_ice(:,:) 564 v_ice_b(:,:) = v_ice(:,:) 565 566 END SUBROUTINE sbc_lim_bef 567 568 SUBROUTINE sbc_lim_diag0 569 !!---------------------------------------------------------------------- 570 !! *** ROUTINE sbc_lim_diag0 *** 571 !! 572 !! ** purpose : set ice-ocean and ice-atm. fluxes to zeros at the beggining 573 !! of the time step 574 !!---------------------------------------------------------------------- 575 sfx (:,:) = 0._wp ; 576 sfx_bri(:,:) = 0._wp ; 577 sfx_sni(:,:) = 0._wp ; sfx_opw(:,:) = 0._wp 578 sfx_bog(:,:) = 0._wp ; sfx_dyn(:,:) = 0._wp 579 sfx_bom(:,:) = 0._wp ; sfx_sum(:,:) = 0._wp 580 sfx_res(:,:) = 0._wp 581 582 wfx_snw(:,:) = 0._wp ; wfx_ice(:,:) = 0._wp 583 wfx_sni(:,:) = 0._wp ; wfx_opw(:,:) = 0._wp 584 wfx_bog(:,:) = 0._wp ; wfx_dyn(:,:) = 0._wp 585 wfx_bom(:,:) = 0._wp ; wfx_sum(:,:) = 0._wp 586 wfx_res(:,:) = 0._wp ; wfx_sub(:,:) = 0._wp 587 wfx_spr(:,:) = 0._wp ; 588 589 hfx_thd(:,:) = 0._wp ; 590 hfx_snw(:,:) = 0._wp ; hfx_opw(:,:) = 0._wp 591 hfx_bog(:,:) = 0._wp ; hfx_dyn(:,:) = 0._wp 592 hfx_bom(:,:) = 0._wp ; hfx_sum(:,:) = 0._wp 593 hfx_res(:,:) = 0._wp ; hfx_sub(:,:) = 0._wp 594 hfx_spr(:,:) = 0._wp ; hfx_dif(:,:) = 0._wp 595 hfx_err(:,:) = 0._wp ; hfx_err_rem(:,:) = 0._wp 596 hfx_err_dif(:,:) = 0._wp ; 597 598 afx_tot(:,:) = 0._wp ; 599 afx_dyn(:,:) = 0._wp ; afx_thd(:,:) = 0._wp 600 601 diag_heat(:,:) = 0._wp ; diag_smvi(:,:) = 0._wp ; 602 diag_vice(:,:) = 0._wp ; diag_vsnw(:,:) = 0._wp ; 603 604 END SUBROUTINE sbc_lim_diag0 605 606 607 FUNCTION fice_cell_ave ( ptab ) 82 608 !!-------------------------------------------------------------------------- 83 609 !! * Compute average over categories, for grid cell (ice covered and free ocean) … … 88 614 89 615 fice_cell_ave (:,:) = 0.0_wp 90 91 616 DO jl = 1, jpl 92 fice_cell_ave (:,:) = fice_cell_ave (:,:) & 93 & + a_i (:,:,jl) * ptab (:,:,jl) 617 fice_cell_ave (:,:) = fice_cell_ave (:,:) + a_i (:,:,jl) * ptab (:,:,jl) 94 618 END DO 95 619 96 620 END FUNCTION fice_cell_ave 97 621 98 FUNCTION fice_ice_ave ( ptab) 622 623 FUNCTION fice_ice_ave ( ptab ) 99 624 !!-------------------------------------------------------------------------- 100 625 !! * Compute average over categories, for ice covered part of grid cell … … 104 629 105 630 fice_ice_ave (:,:) = 0.0_wp 106 WHERE ( at_i (:,:) .GT.0.0_wp ) fice_ice_ave (:,:) = fice_cell_ave ( ptab (:,:,:)) / at_i (:,:)631 WHERE ( at_i (:,:) > 0.0_wp ) fice_ice_ave (:,:) = fice_cell_ave ( ptab (:,:,:)) / at_i (:,:) 107 632 108 633 END FUNCTION fice_ice_ave 109 634 110 !!======================================================================111 112 SUBROUTINE sbc_ice_lim( kt, kblk )113 !!---------------------------------------------------------------------114 !! *** ROUTINE sbc_ice_lim ***115 !!116 !! ** Purpose : update the ocean surface boundary condition via the117 !! Louvain la Neuve Sea Ice Model time stepping118 !!119 !! ** Method : ice model time stepping120 !! - call the ice dynamics routine121 !! - call the ice advection/diffusion routine122 !! - call the ice thermodynamics routine123 !! - call the routine that computes mass and124 !! heat fluxes at the ice/ocean interface125 !! - save the outputs126 !! - save the outputs for restart when necessary127 !!128 !! ** Action : - time evolution of the LIM sea-ice model129 !! - update all sbc variables below sea-ice:130 !! utau, vtau, taum, wndm, qns , qsr, emp , sfx131 !!---------------------------------------------------------------------132 INTEGER, INTENT(in) :: kt ! ocean time step133 INTEGER, INTENT(in) :: kblk ! type of bulk (=3 CLIO, =4 CORE)134 !!135 INTEGER :: jl ! dummy loop index136 REAL(wp) :: zcoef ! local scalar137 REAL(wp), POINTER, DIMENSION(:,:,:) :: zalb_ice_os, zalb_ice_cs ! albedo of the ice under overcast/clear sky138 REAL(wp), POINTER, DIMENSION(:,:,:) :: zalb_ice ! mean albedo of ice (for coupled)139 140 REAL(wp), POINTER, DIMENSION(:,:) :: zalb_ice_all ! Mean albedo over all categories141 REAL(wp), POINTER, DIMENSION(:,:) :: ztem_ice_all ! Mean temperature over all categories142 143 REAL(wp), POINTER, DIMENSION(:,:) :: z_qsr_ice_all ! Mean solar heat flux over all categories144 REAL(wp), POINTER, DIMENSION(:,:) :: z_qns_ice_all ! Mean non solar heat flux over all categories145 REAL(wp), POINTER, DIMENSION(:,:) :: z_qla_ice_all ! Mean latent heat flux over all categories146 REAL(wp), POINTER, DIMENSION(:,:) :: z_dqns_ice_all ! Mean d(qns)/dT over all categories147 REAL(wp), POINTER, DIMENSION(:,:) :: z_dqla_ice_all ! Mean d(qla)/dT over all categories148 !!----------------------------------------------------------------------149 150 !- O.M. : why do we allocate all these arrays even when MOD( kt-1, nn_fsbc ) /= 0 ?????151 152 IF( nn_timing == 1 ) CALL timing_start('sbc_ice_lim')153 154 CALL wrk_alloc( jpi,jpj,jpl, zalb_ice_os, zalb_ice_cs )155 156 #if defined key_coupled157 IF ( ln_cpl .OR. ln_iceflx_ave .OR. ln_iceflx_linear ) CALL wrk_alloc( jpi,jpj,jpl, zalb_ice)158 IF ( ln_iceflx_ave .OR. ln_iceflx_linear ) &159 & CALL wrk_alloc( jpi,jpj, ztem_ice_all, zalb_ice_all, z_qsr_ice_all, z_qns_ice_all, z_qla_ice_all, z_dqns_ice_all, z_dqla_ice_all)160 #endif161 162 IF( kt == nit000 ) THEN163 IF(lwp) WRITE(numout,*)164 IF(lwp) WRITE(numout,*) 'sbc_ice_lim : update ocean surface boudary condition'165 IF(lwp) WRITE(numout,*) '~~~~~~~~~~~ via Louvain la Neuve Ice Model (LIM-3) time stepping'166 !167 CALL ice_init168 !169 IF( ln_nicep ) THEN ! control print at a given point170 jiindx = 177 ; jjindx = 112171 IF(lwp) WRITE(numout,*) ' The debugging point is : jiindx : ',jiindx, ' jjindx : ',jjindx172 ENDIF173 ENDIF174 175 ! !----------------------!176 IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN ! Ice time-step only !177 ! !----------------------!178 ! ! Bulk Formulea !179 ! !----------------!180 !181 u_oce(:,:) = ssu_m(:,:) ! mean surface ocean current at ice velocity point182 v_oce(:,:) = ssv_m(:,:) ! (C-grid dynamics : U- & V-points as the ocean)183 !184 t_bo(:,:) = tfreez( sss_m ) + rt0 ! masked sea surface freezing temperature [Kelvin]185 ! ! (set to rt0 over land)186 CALL albedo_ice( t_su, ht_i, ht_s, zalb_ice_cs, zalb_ice_os ) ! ... ice albedo187 188 DO jl = 1, jpl189 t_su(:,:,jl) = t_su(:,:,jl) + rt0 * ( 1. - tmask(:,:,1) )190 END DO191 192 IF ( ln_cpl ) zalb_ice (:,:,:) = 0.5 * ( zalb_ice_cs (:,:,:) + zalb_ice_os (:,:,:) )193 194 #if defined key_coupled195 IF ( ln_iceflx_ave .OR. ln_iceflx_linear ) THEN196 !197 ! Compute mean albedo and temperature198 zalb_ice_all (:,:) = fice_ice_ave ( zalb_ice (:,:,:) )199 ztem_ice_all (:,:) = fice_ice_ave ( tn_ice (:,:,:) )200 !201 ENDIF202 #endif203 ! Bulk formulea - provides the following fields:204 ! utau_ice, vtau_ice : surface ice stress (U- & V-points) [N/m2]205 ! qsr_ice , qns_ice : solar & non solar heat flux over ice (T-point) [W/m2]206 ! qla_ice : latent heat flux over ice (T-point) [W/m2]207 ! dqns_ice, dqla_ice : non solar & latent heat sensistivity (T-point) [W/m2]208 ! tprecip , sprecip : total & solid precipitation (T-point) [Kg/m2/s]209 ! fr1_i0 , fr2_i0 : 1sr & 2nd fraction of qsr penetration in ice [%]210 !211 SELECT CASE( kblk )212 CASE( 3 ) ! CLIO bulk formulation213 CALL blk_ice_clio( t_su , zalb_ice_cs, zalb_ice_os, &214 & utau_ice , vtau_ice , qns_ice , qsr_ice , &215 & qla_ice , dqns_ice , dqla_ice , &216 & tprecip , sprecip , &217 & fr1_i0 , fr2_i0 , cp_ice_msh, jpl )218 !219 CASE( 4 ) ! CORE bulk formulation220 CALL blk_ice_core( t_su , u_ice , v_ice , zalb_ice_cs, &221 & utau_ice , vtau_ice , qns_ice , qsr_ice , &222 & qla_ice , dqns_ice , dqla_ice , &223 & tprecip , sprecip , &224 & fr1_i0 , fr2_i0 , cp_ice_msh, jpl )225 !226 CASE ( 5 )227 zalb_ice (:,:,:) = 0.5 * ( zalb_ice_cs (:,:,:) + zalb_ice_os (:,:,:) )228 229 CALL sbc_cpl_ice_tau( utau_ice , vtau_ice )230 231 CALL sbc_cpl_ice_flx( p_frld=ato_i, palbi=zalb_ice, psst=sst_m, pist=tn_ice )232 233 ! Latent heat flux is forced to 0 in coupled :234 ! it is included in qns (non-solar heat flux)235 qla_ice (:,:,:) = 0.0e0_wp236 dqla_ice (:,:,:) = 0.0e0_wp237 !238 END SELECT239 240 ! Average over all categories241 #if defined key_coupled242 IF ( ln_iceflx_ave .OR. ln_iceflx_linear ) THEN243 244 z_qns_ice_all (:,:) = fice_ice_ave ( qns_ice (:,:,:) )245 z_qsr_ice_all (:,:) = fice_ice_ave ( qsr_ice (:,:,:) )246 z_dqns_ice_all (:,:) = fice_ice_ave ( dqns_ice (:,:,:) )247 z_qla_ice_all (:,:) = fice_ice_ave ( qla_ice (:,:,:) )248 z_dqla_ice_all (:,:) = fice_ice_ave ( dqla_ice (:,:,:) )249 250 DO jl = 1, jpl251 dqns_ice (:,:,jl) = z_dqns_ice_all (:,:)252 dqla_ice (:,:,jl) = z_dqla_ice_all (:,:)253 END DO254 !255 IF ( ln_iceflx_ave ) THEN256 DO jl = 1, jpl257 qns_ice (:,:,jl) = z_qns_ice_all (:,:)258 qsr_ice (:,:,jl) = z_qsr_ice_all (:,:)259 qla_ice (:,:,jl) = z_qla_ice_all (:,:)260 END DO261 END IF262 !263 IF ( ln_iceflx_linear ) THEN264 DO jl = 1, jpl265 qns_ice (:,:,jl) = z_qns_ice_all(:,:) + z_dqns_ice_all(:,:) * (tn_ice(:,:,jl) - ztem_ice_all(:,:))266 qla_ice (:,:,jl) = z_qla_ice_all(:,:) + z_dqla_ice_all(:,:) * (tn_ice(:,:,jl) - ztem_ice_all(:,:))267 qsr_ice (:,:,jl) = (1.0e0_wp-zalb_ice(:,:,jl)) / (1.0e0_wp-zalb_ice_all(:,:)) * z_qsr_ice_all(:,:)268 END DO269 END IF270 END IF271 #endif272 ! !----------------------!273 ! ! LIM-3 time-stepping !274 ! !----------------------!275 !276 numit = numit + nn_fsbc ! Ice model time step277 !278 ! ! Store previous ice values279 !!gm : remark old_... should becomes ...b as tn versus tb280 old_a_i (:,:,:) = a_i (:,:,:) ! ice area281 old_e_i (:,:,:,:) = e_i (:,:,:,:) ! ice thermal energy282 old_v_i (:,:,:) = v_i (:,:,:) ! ice volume283 old_v_s (:,:,:) = v_s (:,:,:) ! snow volume284 old_e_s (:,:,:,:) = e_s (:,:,:,:) ! snow thermal energy285 old_smv_i(:,:,:) = smv_i(:,:,:) ! salt content286 old_oa_i (:,:,:) = oa_i (:,:,:) ! areal age content287 !288 old_u_ice(:,:) = u_ice(:,:)289 old_v_ice(:,:) = v_ice(:,:)290 ! ! intialisation to zero !!gm is it truly necessary ???291 d_a_i_thd (:,:,:) = 0._wp ; d_a_i_trp (:,:,:) = 0._wp292 d_v_i_thd (:,:,:) = 0._wp ; d_v_i_trp (:,:,:) = 0._wp293 d_e_i_thd (:,:,:,:) = 0._wp ; d_e_i_trp (:,:,:,:) = 0._wp294 d_v_s_thd (:,:,:) = 0._wp ; d_v_s_trp (:,:,:) = 0._wp295 d_e_s_thd (:,:,:,:) = 0._wp ; d_e_s_trp (:,:,:,:) = 0._wp296 d_smv_i_thd(:,:,:) = 0._wp ; d_smv_i_trp(:,:,:) = 0._wp297 d_oa_i_thd (:,:,:) = 0._wp ; d_oa_i_trp (:,:,:) = 0._wp298 !299 d_u_ice_dyn(:,:) = 0._wp300 d_v_ice_dyn(:,:) = 0._wp301 !302 sfx (:,:) = 0._wp ; sfx_thd (:,:) = 0._wp303 sfx_bri(:,:) = 0._wp ; sfx_mec (:,:) = 0._wp ; sfx_res (:,:) = 0._wp304 fhbri (:,:) = 0._wp ; fheat_mec(:,:) = 0._wp ; fheat_res(:,:) = 0._wp305 fhmec (:,:) = 0._wp ;306 fmmec (:,:) = 0._wp307 fmmflx (:,:) = 0._wp308 focea2D(:,:) = 0._wp309 fsup2D (:,:) = 0._wp310 311 ! used in limthd.F90312 rdvosif(:,:) = 0._wp ! variation of ice volume at surface313 rdvobif(:,:) = 0._wp ! variation of ice volume at bottom314 fdvolif(:,:) = 0._wp ! total variation of ice volume315 rdvonif(:,:) = 0._wp ! lateral variation of ice volume316 fstric (:,:) = 0._wp ! part of solar radiation transmitted through the ice317 ffltbif(:,:) = 0._wp ! linked with fstric318 qfvbq (:,:) = 0._wp ! linked with fstric319 rdm_snw(:,:) = 0._wp ! variation of snow mass per unit area320 rdm_ice(:,:) = 0._wp ! variation of ice mass per unit area321 hicifp (:,:) = 0._wp ! daily thermodynamic ice production.322 !323 diag_sni_gr(:,:) = 0._wp ; diag_lat_gr(:,:) = 0._wp324 diag_bot_gr(:,:) = 0._wp ; diag_dyn_gr(:,:) = 0._wp325 diag_bot_me(:,:) = 0._wp ; diag_sur_me(:,:) = 0._wp326 diag_res_pr(:,:) = 0._wp ; diag_trp_vi(:,:) = 0._wp327 ! dynamical invariants328 delta_i(:,:) = 0._wp ; divu_i(:,:) = 0._wp ; shear_i(:,:) = 0._wp329 330 CALL lim_rst_opn( kt ) ! Open Ice restart file331 !332 IF( ln_nicep ) CALL lim_prt_state( kt, jiindx, jjindx, 1, ' - Beginning the time step - ' ) ! control print333 ! ----------------------------------------------334 ! ice dynamics and transport (except in 1D case)335 ! ----------------------------------------------336 IF( .NOT. lk_c1d ) THEN337 CALL lim_dyn( kt ) ! Ice dynamics ( rheology/dynamics )338 CALL lim_trp( kt ) ! Ice transport ( Advection/diffusion )339 CALL lim_var_glo2eqv ! equivalent variables, requested for rafting340 IF( ln_nicep ) CALL lim_prt_state( kt, jiindx, jjindx,-1, ' - ice dyn & trp - ' ) ! control print341 CALL lim_itd_me ! Mechanical redistribution ! (ridging/rafting)342 CALL lim_var_agg( 1 )343 #if defined key_bdy344 ! bdy ice thermo345 CALL lim_var_glo2eqv ! equivalent variables346 CALL bdy_ice_lim( kt )347 CALL lim_itd_me_zapsmall348 CALL lim_var_agg(1)349 IF( ln_nicep ) CALL lim_prt_state( kt, jiindx, jjindx, 1, ' - ice thermo bdy - ' ) ! control print350 #endif351 CALL lim_update1352 ENDIF353 ! !- Change old values for new values354 old_u_ice(:,:) = u_ice (:,:)355 old_v_ice(:,:) = v_ice (:,:)356 old_a_i(:,:,:) = a_i (:,:,:)357 old_v_s(:,:,:) = v_s (:,:,:)358 old_v_i(:,:,:) = v_i (:,:,:)359 old_e_s(:,:,:,:) = e_s (:,:,:,:)360 old_e_i(:,:,:,:) = e_i (:,:,:,:)361 old_oa_i(:,:,:) = oa_i(:,:,:)362 old_smv_i(:,:,:) = smv_i (:,:,:)363 364 ! ----------------------------------------------365 ! ice thermodynamic366 ! ----------------------------------------------367 CALL lim_var_glo2eqv ! equivalent variables368 CALL lim_var_agg(1) ! aggregate ice categories369 ! previous lead fraction and ice volume for flux calculations370 pfrld(:,:) = 1._wp - at_i(:,:)371 phicif(:,:) = vt_i(:,:)372 !373 CALL lim_var_bv ! bulk brine volume (diag)374 CALL lim_thd( kt ) ! Ice thermodynamics375 zcoef = rdt_ice /rday ! Ice natural aging376 oa_i(:,:,:) = oa_i(:,:,:) + a_i(:,:,:) * zcoef377 CALL lim_var_glo2eqv ! this CALL is maybe not necessary (Martin)378 IF( ln_nicep ) CALL lim_prt_state( kt, jiindx, jjindx, 1, ' - ice thermodyn. - ' ) ! control print379 CALL lim_itd_th( kt ) ! Remap ice categories, lateral accretion !380 CALL lim_var_agg( 1 ) ! requested by limupdate381 CALL lim_update2 ! Global variables update382 383 CALL lim_var_glo2eqv ! equivalent variables (outputs)384 CALL lim_var_agg(2) ! aggregate ice thickness categories385 IF( ln_nicep ) CALL lim_prt_state( kt, jiindx, jjindx, 2, ' - Final state - ' ) ! control print386 !387 CALL lim_sbc_flx( kt ) ! Update surface ocean mass, heat and salt fluxes388 !389 IF( ln_nicep ) CALL lim_prt_state( kt, jiindx, jjindx, 3, ' - Final state lim_sbc - ' ) ! control print390 !391 ! ! Diagnostics and outputs392 IF (ln_limdiaout) CALL lim_diahsb393 !clem # if ! defined key_iomput394 CALL lim_wri( 1 ) ! Ice outputs395 !clem # endif396 IF( kt == nit000 .AND. ln_rstart ) &397 & CALL iom_close( numrir ) ! clem: close input ice restart file398 !399 IF( lrst_ice ) CALL lim_rst_write( kt ) ! Ice restart file400 CALL lim_var_glo2eqv ! ???401 !402 IF( ln_nicep ) CALL lim_ctl( kt ) ! alerts in case of model crash403 !404 ENDIF ! End sea-ice time step only405 406 ! !--------------------------!407 ! ! at all ocean time step !408 ! !--------------------------!409 !410 ! ! Update surface ocean stresses (only in ice-dynamic case)411 ! ! otherwise the atm.-ocean stresses are used everywhere412 IF( ln_limdyn ) CALL lim_sbc_tau( kt, ub(:,:,1), vb(:,:,1) ) ! using before instantaneous surf. currents413 414 !!gm remark, the ocean-ice stress is not saved in ice diag call above ..... find a solution!!!415 !416 CALL wrk_dealloc( jpi,jpj,jpl, zalb_ice_os, zalb_ice_cs )417 418 #if defined key_coupled419 IF ( ln_cpl .OR. ln_iceflx_ave .OR. ln_iceflx_linear ) CALL wrk_dealloc( jpi,jpj,jpl, zalb_ice)420 IF ( ln_iceflx_ave .OR. ln_iceflx_linear ) &421 & CALL wrk_dealloc( jpi,jpj, ztem_ice_all, zalb_ice_all, z_qsr_ice_all, z_qns_ice_all, z_qla_ice_all, z_dqns_ice_all, z_dqla_ice_all)422 #endif423 !424 IF( nn_timing == 1 ) CALL timing_stop('sbc_ice_lim')425 !426 END SUBROUTINE sbc_ice_lim427 428 429 SUBROUTINE lim_ctl( kt )430 !!-----------------------------------------------------------------------431 !! *** ROUTINE lim_ctl ***432 !!433 !! ** Purpose : Alerts in case of model crash434 !!-------------------------------------------------------------------435 INTEGER, INTENT(in) :: kt ! ocean time step436 INTEGER :: ji, jj, jk, jl ! dummy loop indices437 INTEGER :: inb_altests ! number of alert tests (max 20)438 INTEGER :: ialert_id ! number of the current alert439 REAL(wp) :: ztmelts ! ice layer melting point440 CHARACTER (len=30), DIMENSION(20) :: cl_alname ! name of alert441 INTEGER , DIMENSION(20) :: inb_alp ! number of alerts positive442 !!-------------------------------------------------------------------443 444 inb_altests = 10445 inb_alp(:) = 0446 447 ! Alert if incompatible volume and concentration448 ialert_id = 2 ! reference number of this alert449 cl_alname(ialert_id) = ' Incompat vol and con ' ! name of the alert450 451 DO jl = 1, jpl452 DO jj = 1, jpj453 DO ji = 1, jpi454 IF( v_i(ji,jj,jl) /= 0._wp .AND. a_i(ji,jj,jl) == 0._wp ) THEN455 !WRITE(numout,*) ' ALERTE 2 : Incompatible volume and concentration '456 !WRITE(numout,*) ' at_i ', at_i(ji,jj)457 !WRITE(numout,*) ' Point - category', ji, jj, jl458 !WRITE(numout,*) ' a_i *** a_i_old ', a_i (ji,jj,jl), old_a_i (ji,jj,jl)459 !WRITE(numout,*) ' v_i *** v_i_old ', v_i (ji,jj,jl), old_v_i (ji,jj,jl)460 !WRITE(numout,*) ' d_a_i_thd/trp ', d_a_i_thd(ji,jj,jl), d_a_i_trp(ji,jj,jl)461 !WRITE(numout,*) ' d_v_i_thd/trp ', d_v_i_thd(ji,jj,jl), d_v_i_trp(ji,jj,jl)462 inb_alp(ialert_id) = inb_alp(ialert_id) + 1463 ENDIF464 END DO465 END DO466 END DO467 468 ! Alerte if very thick ice469 ialert_id = 3 ! reference number of this alert470 cl_alname(ialert_id) = ' Very thick ice ' ! name of the alert471 jl = jpl472 DO jj = 1, jpj473 DO ji = 1, jpi474 IF( ht_i(ji,jj,jl) > 50._wp ) THEN475 !CALL lim_prt_state( kt, ji, jj, 2, ' ALERTE 3 : Very thick ice ' )476 inb_alp(ialert_id) = inb_alp(ialert_id) + 1477 ENDIF478 END DO479 END DO480 481 ! Alert if very fast ice482 ialert_id = 4 ! reference number of this alert483 cl_alname(ialert_id) = ' Very fast ice ' ! name of the alert484 DO jj = 1, jpj485 DO ji = 1, jpi486 IF( MAX( ABS( u_ice(ji,jj) ), ABS( v_ice(ji,jj) ) ) > 1.5 .AND. &487 & at_i(ji,jj) > 0._wp ) THEN488 !CALL lim_prt_state( kt, ji, jj, 1, ' ALERTE 4 : Very fast ice ' )489 !WRITE(numout,*) ' ice strength : ', strength(ji,jj)490 !WRITE(numout,*) ' oceanic stress utau : ', utau(ji,jj)491 !WRITE(numout,*) ' oceanic stress vtau : ', vtau(ji,jj)492 !WRITE(numout,*) ' sea-ice stress utau_ice : ', utau_ice(ji,jj)493 !WRITE(numout,*) ' sea-ice stress vtau_ice : ', vtau_ice(ji,jj)494 !WRITE(numout,*) ' oceanic speed u : ', u_oce(ji,jj)495 !WRITE(numout,*) ' oceanic speed v : ', v_oce(ji,jj)496 !WRITE(numout,*) ' sst : ', sst_m(ji,jj)497 !WRITE(numout,*) ' sss : ', sss_m(ji,jj)498 !WRITE(numout,*)499 inb_alp(ialert_id) = inb_alp(ialert_id) + 1500 ENDIF501 END DO502 END DO503 504 ! Alert if there is ice on continents505 ialert_id = 6 ! reference number of this alert506 cl_alname(ialert_id) = ' Ice on continents ' ! name of the alert507 DO jj = 1, jpj508 DO ji = 1, jpi509 IF( tms(ji,jj) <= 0._wp .AND. at_i(ji,jj) > 0._wp ) THEN510 !CALL lim_prt_state( kt, ji, jj, 1, ' ALERTE 6 : Ice on continents ' )511 !WRITE(numout,*) ' masks s, u, v : ', tms(ji,jj), tmu(ji,jj), tmv(ji,jj)512 !WRITE(numout,*) ' sst : ', sst_m(ji,jj)513 !WRITE(numout,*) ' sss : ', sss_m(ji,jj)514 !WRITE(numout,*) ' at_i(ji,jj) : ', at_i(ji,jj)515 !WRITE(numout,*) ' v_ice(ji,jj) : ', v_ice(ji,jj)516 !WRITE(numout,*) ' v_ice(ji,jj-1) : ', v_ice(ji,jj-1)517 !WRITE(numout,*) ' u_ice(ji-1,jj) : ', u_ice(ji-1,jj)518 !WRITE(numout,*) ' u_ice(ji,jj) : ', v_ice(ji,jj)519 !520 inb_alp(ialert_id) = inb_alp(ialert_id) + 1521 ENDIF522 END DO523 END DO524 525 !526 ! ! Alert if very fresh ice527 ialert_id = 7 ! reference number of this alert528 cl_alname(ialert_id) = ' Very fresh ice ' ! name of the alert529 DO jl = 1, jpl530 DO jj = 1, jpj531 DO ji = 1, jpi532 IF( sm_i(ji,jj,jl) < 0.1 .AND. a_i(ji,jj,jl) > 0._wp ) THEN533 ! CALL lim_prt_state(kt,ji,jj,1, ' ALERTE 7 : Very fresh ice ' )534 ! WRITE(numout,*) ' sst : ', sst_m(ji,jj)535 ! WRITE(numout,*) ' sss : ', sss_m(ji,jj)536 ! WRITE(numout,*) ' s_i_newice : ', s_i_newice(ji,jj,1:jpl)537 ! WRITE(numout,*)538 inb_alp(ialert_id) = inb_alp(ialert_id) + 1539 ENDIF540 END DO541 END DO542 END DO543 !544 545 ! ! Alert if too old ice546 ialert_id = 9 ! reference number of this alert547 cl_alname(ialert_id) = ' Very old ice ' ! name of the alert548 DO jl = 1, jpl549 DO jj = 1, jpj550 DO ji = 1, jpi551 IF ( ( ( ABS( o_i(ji,jj,jl) ) > rdt_ice ) .OR. &552 ( ABS( o_i(ji,jj,jl) ) < 0._wp) ) .AND. &553 ( a_i(ji,jj,jl) > 0._wp ) ) THEN554 !CALL lim_prt_state( kt, ji, jj, 1, ' ALERTE 9 : Wrong ice age ')555 inb_alp(ialert_id) = inb_alp(ialert_id) + 1556 ENDIF557 END DO558 END DO559 END DO560 561 ! Alert on salt flux562 ialert_id = 5 ! reference number of this alert563 cl_alname(ialert_id) = ' High salt flux ' ! name of the alert564 DO jj = 1, jpj565 DO ji = 1, jpi566 IF( ABS( sfx (ji,jj) ) .GT. 1.0e-2 ) THEN ! = 1 psu/day for 1m ocean depth567 !CALL lim_prt_state( kt, ji, jj, 3, ' ALERTE 5 : High salt flux ' )568 !DO jl = 1, jpl569 !WRITE(numout,*) ' Category no: ', jl570 !WRITE(numout,*) ' a_i : ', a_i (ji,jj,jl) , ' old_a_i : ', old_a_i (ji,jj,jl)571 !WRITE(numout,*) ' d_a_i_trp : ', d_a_i_trp(ji,jj,jl) , ' d_a_i_thd : ', d_a_i_thd(ji,jj,jl)572 !WRITE(numout,*) ' v_i : ', v_i (ji,jj,jl) , ' old_v_i : ', old_v_i (ji,jj,jl)573 !WRITE(numout,*) ' d_v_i_trp : ', d_v_i_trp(ji,jj,jl) , ' d_v_i_thd : ', d_v_i_thd(ji,jj,jl)574 !WRITE(numout,*) ' '575 !END DO576 inb_alp(ialert_id) = inb_alp(ialert_id) + 1577 ENDIF578 END DO579 END DO580 581 ! Alert if qns very big582 ialert_id = 8 ! reference number of this alert583 cl_alname(ialert_id) = ' fnsolar very big ' ! name of the alert584 DO jj = 1, jpj585 DO ji = 1, jpi586 IF( ABS( qns(ji,jj) ) > 1500._wp .AND. at_i(ji,jj) > 0._wp ) THEN587 !588 !WRITE(numout,*) ' ALERTE 8 : Very high non-solar heat flux'589 !WRITE(numout,*) ' ji, jj : ', ji, jj590 !WRITE(numout,*) ' qns : ', qns(ji,jj)591 !WRITE(numout,*) ' sst : ', sst_m(ji,jj)592 !WRITE(numout,*) ' sss : ', sss_m(ji,jj)593 !WRITE(numout,*) ' qcmif : ', qcmif(ji,jj)594 !WRITE(numout,*) ' qldif : ', qldif(ji,jj)595 !WRITE(numout,*) ' qcmif : ', qcmif(ji,jj) / rdt_ice596 !WRITE(numout,*) ' qldif : ', qldif(ji,jj) / rdt_ice597 !WRITE(numout,*) ' qfvbq : ', qfvbq(ji,jj)598 !WRITE(numout,*) ' qdtcn : ', qdtcn(ji,jj)599 !WRITE(numout,*) ' qfvbq / dt: ', qfvbq(ji,jj) / rdt_ice600 !WRITE(numout,*) ' qdtcn / dt: ', qdtcn(ji,jj) / rdt_ice601 !WRITE(numout,*) ' fdtcn : ', fdtcn(ji,jj)602 !WRITE(numout,*) ' fhmec : ', fhmec(ji,jj)603 !WRITE(numout,*) ' fheat_mec : ', fheat_mec(ji,jj)604 !WRITE(numout,*) ' fheat_res : ', fheat_res(ji,jj)605 !WRITE(numout,*) ' fhbri : ', fhbri(ji,jj)606 !607 !CALL lim_prt_state( kt, ji, jj, 2, ' ')608 inb_alp(ialert_id) = inb_alp(ialert_id) + 1609 !610 ENDIF611 END DO612 END DO613 !+++++614 615 ! Alert if very warm ice616 ialert_id = 10 ! reference number of this alert617 cl_alname(ialert_id) = ' Very warm ice ' ! name of the alert618 inb_alp(ialert_id) = 0619 DO jl = 1, jpl620 DO jk = 1, nlay_i621 DO jj = 1, jpj622 DO ji = 1, jpi623 ztmelts = -tmut * s_i(ji,jj,jk,jl) + rtt624 IF( t_i(ji,jj,jk,jl) >= ztmelts .AND. v_i(ji,jj,jl) > 1.e-10 &625 & .AND. a_i(ji,jj,jl) > 0._wp ) THEN626 !WRITE(numout,*) ' ALERTE 10 : Very warm ice'627 !WRITE(numout,*) ' ji, jj, jk, jl : ', ji, jj, jk, jl628 !WRITE(numout,*) ' t_i : ', t_i(ji,jj,jk,jl)629 !WRITE(numout,*) ' e_i : ', e_i(ji,jj,jk,jl)630 !WRITE(numout,*) ' s_i : ', s_i(ji,jj,jk,jl)631 !WRITE(numout,*) ' ztmelts : ', ztmelts632 inb_alp(ialert_id) = inb_alp(ialert_id) + 1633 ENDIF634 END DO635 END DO636 END DO637 END DO638 639 ! sum of the alerts on all processors640 IF( lk_mpp ) THEN641 DO ialert_id = 1, inb_altests642 CALL mpp_sum(inb_alp(ialert_id))643 END DO644 ENDIF645 646 ! print alerts647 IF( lwp ) THEN648 ialert_id = 1 ! reference number of this alert649 cl_alname(ialert_id) = ' NO alerte 1 ' ! name of the alert650 WRITE(numout,*) ' time step ',kt651 WRITE(numout,*) ' All alerts at the end of ice model '652 DO ialert_id = 1, inb_altests653 WRITE(numout,*) ialert_id, cl_alname(ialert_id)//' : ', inb_alp(ialert_id), ' times ! '654 END DO655 ENDIF656 !657 END SUBROUTINE lim_ctl658 659 660 SUBROUTINE lim_prt_state( kt, ki, kj, kn, cd1 )661 !!-----------------------------------------------------------------------662 !! *** ROUTINE lim_prt_state ***663 !!664 !! ** Purpose : Writes global ice state on the (i,j) point665 !! in ocean.ouput666 !! 3 possibilities exist667 !! n = 1/-1 -> simple ice state (plus Mechanical Check if -1)668 !! n = 2 -> exhaustive state669 !! n = 3 -> ice/ocean salt fluxes670 !!671 !! ** input : point coordinates (i,j)672 !! n : number of the option673 !!-------------------------------------------------------------------674 INTEGER , INTENT(in) :: kt ! ocean time step675 INTEGER , INTENT(in) :: ki, kj, kn ! ocean gridpoint indices676 CHARACTER(len=*), INTENT(in) :: cd1 !677 !!678 INTEGER :: jl, ji, jj679 !!-------------------------------------------------------------------680 681 DO ji = mi0(ki), mi1(ki)682 DO jj = mj0(kj), mj1(kj)683 684 WRITE(numout,*) ' time step ',kt,' ',cd1 ! print title685 686 !----------------687 ! Simple state688 !----------------689 690 IF ( kn == 1 .OR. kn == -1 ) THEN691 WRITE(numout,*) ' lim_prt_state - Point : ',ji,jj692 WRITE(numout,*) ' ~~~~~~~~~~~~~~ '693 WRITE(numout,*) ' Simple state '694 WRITE(numout,*) ' masks s,u,v : ', tms(ji,jj), tmu(ji,jj), tmv(ji,jj)695 WRITE(numout,*) ' lat - long : ', gphit(ji,jj), glamt(ji,jj)696 WRITE(numout,*) ' Time step : ', numit697 WRITE(numout,*) ' - Ice drift '698 WRITE(numout,*) ' ~~~~~~~~~~~ '699 WRITE(numout,*) ' u_ice(i-1,j) : ', u_ice(ji-1,jj)700 WRITE(numout,*) ' u_ice(i ,j) : ', u_ice(ji,jj)701 WRITE(numout,*) ' v_ice(i ,j-1): ', v_ice(ji,jj-1)702 WRITE(numout,*) ' v_ice(i ,j) : ', v_ice(ji,jj)703 WRITE(numout,*) ' strength : ', strength(ji,jj)704 WRITE(numout,*)705 WRITE(numout,*) ' - Cell values '706 WRITE(numout,*) ' ~~~~~~~~~~~ '707 WRITE(numout,*) ' cell area : ', area(ji,jj)708 WRITE(numout,*) ' at_i : ', at_i(ji,jj)709 WRITE(numout,*) ' vt_i : ', vt_i(ji,jj)710 WRITE(numout,*) ' vt_s : ', vt_s(ji,jj)711 DO jl = 1, jpl712 WRITE(numout,*) ' - Category (', jl,')'713 WRITE(numout,*) ' a_i : ', a_i(ji,jj,jl)714 WRITE(numout,*) ' ht_i : ', ht_i(ji,jj,jl)715 WRITE(numout,*) ' ht_s : ', ht_s(ji,jj,jl)716 WRITE(numout,*) ' v_i : ', v_i(ji,jj,jl)717 WRITE(numout,*) ' v_s : ', v_s(ji,jj,jl)718 WRITE(numout,*) ' e_s : ', e_s(ji,jj,1,jl)/1.0e9719 WRITE(numout,*) ' e_i : ', e_i(ji,jj,1:nlay_i,jl)/1.0e9720 WRITE(numout,*) ' t_su : ', t_su(ji,jj,jl)721 WRITE(numout,*) ' t_snow : ', t_s(ji,jj,1,jl)722 WRITE(numout,*) ' t_i : ', t_i(ji,jj,1:nlay_i,jl)723 WRITE(numout,*) ' sm_i : ', sm_i(ji,jj,jl)724 WRITE(numout,*) ' smv_i : ', smv_i(ji,jj,jl)725 WRITE(numout,*)726 END DO727 ENDIF728 IF( kn == -1 ) THEN729 WRITE(numout,*) ' Mechanical Check ************** '730 WRITE(numout,*) ' Check what means ice divergence '731 WRITE(numout,*) ' Total ice concentration ', at_i (ji,jj)732 WRITE(numout,*) ' Total lead fraction ', ato_i(ji,jj)733 WRITE(numout,*) ' Sum of both ', ato_i(ji,jj) + at_i(ji,jj)734 WRITE(numout,*) ' Sum of both minus 1 ', ato_i(ji,jj) + at_i(ji,jj) - 1.00735 ENDIF736 737 738 !--------------------739 ! Exhaustive state740 !--------------------741 742 IF ( kn .EQ. 2 ) THEN743 WRITE(numout,*) ' lim_prt_state - Point : ',ji,jj744 WRITE(numout,*) ' ~~~~~~~~~~~~~~ '745 WRITE(numout,*) ' Exhaustive state '746 WRITE(numout,*) ' lat - long ', gphit(ji,jj), glamt(ji,jj)747 WRITE(numout,*) ' Time step ', numit748 WRITE(numout,*)749 WRITE(numout,*) ' - Cell values '750 WRITE(numout,*) ' ~~~~~~~~~~~ '751 WRITE(numout,*) ' cell area : ', area(ji,jj)752 WRITE(numout,*) ' at_i : ', at_i(ji,jj)753 WRITE(numout,*) ' vt_i : ', vt_i(ji,jj)754 WRITE(numout,*) ' vt_s : ', vt_s(ji,jj)755 WRITE(numout,*) ' u_ice(i-1,j) : ', u_ice(ji-1,jj)756 WRITE(numout,*) ' u_ice(i ,j) : ', u_ice(ji,jj)757 WRITE(numout,*) ' v_ice(i ,j-1): ', v_ice(ji,jj-1)758 WRITE(numout,*) ' v_ice(i ,j) : ', v_ice(ji,jj)759 WRITE(numout,*) ' strength : ', strength(ji,jj)760 WRITE(numout,*) ' d_u_ice_dyn : ', d_u_ice_dyn(ji,jj), ' d_v_ice_dyn : ', d_v_ice_dyn(ji,jj)761 WRITE(numout,*) ' old_u_ice : ', old_u_ice(ji,jj) , ' old_v_ice : ', old_v_ice(ji,jj)762 WRITE(numout,*)763 764 DO jl = 1, jpl765 WRITE(numout,*) ' - Category (',jl,')'766 WRITE(numout,*) ' ~~~~~~~~ '767 WRITE(numout,*) ' ht_i : ', ht_i(ji,jj,jl) , ' ht_s : ', ht_s(ji,jj,jl)768 WRITE(numout,*) ' t_i : ', t_i(ji,jj,1:nlay_i,jl)769 WRITE(numout,*) ' t_su : ', t_su(ji,jj,jl) , ' t_s : ', t_s(ji,jj,1,jl)770 WRITE(numout,*) ' sm_i : ', sm_i(ji,jj,jl) , ' o_i : ', o_i(ji,jj,jl)771 WRITE(numout,*) ' a_i : ', a_i(ji,jj,jl) , ' old_a_i : ', old_a_i(ji,jj,jl)772 WRITE(numout,*) ' d_a_i_trp : ', d_a_i_trp(ji,jj,jl) , ' d_a_i_thd : ', d_a_i_thd(ji,jj,jl)773 WRITE(numout,*) ' v_i : ', v_i(ji,jj,jl) , ' old_v_i : ', old_v_i(ji,jj,jl)774 WRITE(numout,*) ' d_v_i_trp : ', d_v_i_trp(ji,jj,jl) , ' d_v_i_thd : ', d_v_i_thd(ji,jj,jl)775 WRITE(numout,*) ' v_s : ', v_s(ji,jj,jl) , ' old_v_s : ', old_v_s(ji,jj,jl)776 WRITE(numout,*) ' d_v_s_trp : ', d_v_s_trp(ji,jj,jl) , ' d_v_s_thd : ', d_v_s_thd(ji,jj,jl)777 WRITE(numout,*) ' e_i1 : ', e_i(ji,jj,1,jl)/1.0e9 , ' old_ei1 : ', old_e_i(ji,jj,1,jl)/1.0e9778 WRITE(numout,*) ' de_i1_trp : ', d_e_i_trp(ji,jj,1,jl)/1.0e9, ' de_i1_thd : ', d_e_i_thd(ji,jj,1,jl)/1.0e9779 WRITE(numout,*) ' e_i2 : ', e_i(ji,jj,2,jl)/1.0e9 , ' old_ei2 : ', old_e_i(ji,jj,2,jl)/1.0e9780 WRITE(numout,*) ' de_i2_trp : ', d_e_i_trp(ji,jj,2,jl)/1.0e9, ' de_i2_thd : ', d_e_i_thd(ji,jj,2,jl)/1.0e9781 WRITE(numout,*) ' e_snow : ', e_s(ji,jj,1,jl) , ' old_e_snow : ', old_e_s(ji,jj,1,jl)782 WRITE(numout,*) ' d_e_s_trp : ', d_e_s_trp(ji,jj,1,jl) , ' d_e_s_thd : ', d_e_s_thd(ji,jj,1,jl)783 WRITE(numout,*) ' smv_i : ', smv_i(ji,jj,jl) , ' old_smv_i : ', old_smv_i(ji,jj,jl)784 WRITE(numout,*) ' d_smv_i_trp: ', d_smv_i_trp(ji,jj,jl) , ' d_smv_i_thd: ', d_smv_i_thd(ji,jj,jl)785 WRITE(numout,*) ' oa_i : ', oa_i(ji,jj,jl) , ' old_oa_i : ', old_oa_i(ji,jj,jl)786 WRITE(numout,*) ' d_oa_i_trp : ', d_oa_i_trp(ji,jj,jl) , ' d_oa_i_thd : ', d_oa_i_thd(ji,jj,jl)787 END DO !jl788 789 WRITE(numout,*)790 WRITE(numout,*) ' - Heat / FW fluxes '791 WRITE(numout,*) ' ~~~~~~~~~~~~~~~~ '792 WRITE(numout,*) ' emp : ', emp (ji,jj)793 WRITE(numout,*) ' sfx : ', sfx (ji,jj)794 WRITE(numout,*) ' sfx_thd : ', sfx_thd(ji,jj)795 WRITE(numout,*) ' sfx_bri : ', sfx_bri (ji,jj)796 WRITE(numout,*) ' sfx_mec : ', sfx_mec (ji,jj)797 WRITE(numout,*) ' sfx_res : ', sfx_res(ji,jj)798 WRITE(numout,*) ' fmmec : ', fmmec (ji,jj)799 WRITE(numout,*) ' fhmec : ', fhmec (ji,jj)800 WRITE(numout,*) ' fhbri : ', fhbri (ji,jj)801 WRITE(numout,*) ' fheat_mec : ', fheat_mec(ji,jj)802 WRITE(numout,*)803 WRITE(numout,*) ' sst : ', sst_m(ji,jj)804 WRITE(numout,*) ' sss : ', sss_m(ji,jj)805 WRITE(numout,*)806 WRITE(numout,*) ' - Stresses '807 WRITE(numout,*) ' ~~~~~~~~ '808 WRITE(numout,*) ' utau_ice : ', utau_ice(ji,jj)809 WRITE(numout,*) ' vtau_ice : ', vtau_ice(ji,jj)810 WRITE(numout,*) ' utau : ', utau (ji,jj)811 WRITE(numout,*) ' vtau : ', vtau (ji,jj)812 WRITE(numout,*) ' oc. vel. u : ', u_oce (ji,jj)813 WRITE(numout,*) ' oc. vel. v : ', v_oce (ji,jj)814 ENDIF815 816 !---------------------817 ! Salt / heat fluxes818 !---------------------819 820 IF ( kn .EQ. 3 ) THEN821 WRITE(numout,*) ' lim_prt_state - Point : ',ji,jj822 WRITE(numout,*) ' ~~~~~~~~~~~~~~ '823 WRITE(numout,*) ' - Salt / Heat Fluxes '824 WRITE(numout,*) ' ~~~~~~~~~~~~~~~~ '825 WRITE(numout,*) ' lat - long ', gphit(ji,jj), glamt(ji,jj)826 WRITE(numout,*) ' Time step ', numit827 WRITE(numout,*)828 WRITE(numout,*) ' - Heat fluxes at bottom interface ***'829 WRITE(numout,*) ' qsr : ', qsr(ji,jj)830 WRITE(numout,*) ' qns : ', qns(ji,jj)831 WRITE(numout,*) ' fdtcn : ', fdtcn(ji,jj)832 WRITE(numout,*) ' qcmif : ', qcmif(ji,jj) * r1_rdtice833 WRITE(numout,*) ' qldif : ', qldif(ji,jj) * r1_rdtice834 WRITE(numout,*)835 WRITE(numout,*) ' - Salt fluxes at bottom interface ***'836 WRITE(numout,*) ' emp : ', emp (ji,jj)837 WRITE(numout,*) ' sfx_bri : ', sfx_bri(ji,jj)838 WRITE(numout,*) ' sfx : ', sfx (ji,jj)839 WRITE(numout,*) ' sfx_res : ', sfx_res(ji,jj)840 WRITE(numout,*) ' sfx_mec : ', sfx_mec(ji,jj)841 WRITE(numout,*) ' - Heat fluxes at bottom interface ***'842 WRITE(numout,*) ' fheat_res : ', fheat_res(ji,jj)843 WRITE(numout,*)844 WRITE(numout,*) ' - Momentum fluxes '845 WRITE(numout,*) ' utau : ', utau(ji,jj)846 WRITE(numout,*) ' vtau : ', vtau(ji,jj)847 ENDIF848 WRITE(numout,*) ' '849 !850 END DO851 END DO852 853 END SUBROUTINE lim_prt_state854 635 855 636 #else … … 861 642 WRITE(*,*) 'sbc_ice_lim: You should not have seen this print! error?', kt, kblk 862 643 END SUBROUTINE sbc_ice_lim 644 SUBROUTINE sbc_lim_init ! Dummy routine 645 END SUBROUTINE sbc_lim_init 863 646 #endif 864 647
Note: See TracChangeset
for help on using the changeset viewer.