Changeset 1239
- Timestamp:
- 2008-12-31T10:35:35+01:00 (15 years ago)
- Location:
- trunk/NEMO
- Files:
-
- 1 added
- 9 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/NEMO/OPA_SRC/DYN/dynzdf.F90
r1152 r1239 106 106 !!---------------------------------------------------------------------- 107 107 USE zdftke 108 USE zdftke2 108 109 USE zdfkpp 109 110 !!---------------------------------------------------------------------- … … 115 116 116 117 ! Force implicit schemes 117 IF( lk_zdftke .OR. lk_zdf kpp )nzdf = 1 ! TKE or KPP physics118 IF( ln_dynldf_iso ) nzdf = 1 ! iso-neutral lateral physics119 IF( ln_dynldf_hor .AND. ln_sco ) nzdf = 1 ! horizontal lateral physics in s-coordinate118 IF( lk_zdftke .OR. lk_zdftke2 .OR. lk_zdfkpp ) nzdf = 1 ! TKE or KPP physics 119 IF( ln_dynldf_iso ) nzdf = 1 ! iso-neutral lateral physics 120 IF( ln_dynldf_hor .AND. ln_sco ) nzdf = 1 ! horizontal lateral physics in s-coordinate 120 121 121 122 IF( lk_esopa ) nzdf = -1 ! Esopa key: All schemes used -
trunk/NEMO/OPA_SRC/IOM/in_out_manager.F90
r1229 r1239 43 43 !!---------------------------------------------------------------------- 44 44 INTEGER :: nitrst !: time step at which restart file should be written 45 #if defined key_zdftke2 46 INTEGER :: nitrst_tke2 !: time step at which restart file should be written 47 #endif 45 48 !!---------------------------------------------------------------------- 46 49 !! output monitoring -
trunk/NEMO/OPA_SRC/SBC/sbcice_lim.F90
r1220 r1239 233 233 ! 234 234 ! ! Diagnostics and outputs 235 ! ! Ice Diagnostics 235 #if defined key_zdftke2 236 IF( kt .NE. nitend+1 ) THEN 237 IF( ( MOD( kt+nn_fsbc-1, ninfo ) == 0 .OR. ntmoy == 1 ) .AND. .NOT. lk_mpp ) & 238 & CALL lim_dia 239 CALL lim_wri( 1 ) ! Ice outputs 240 IF( lrst_ice ) CALL lim_rst_write( kt ) ! Ice restart file 241 CALL lim_var_glo2eqv ! ??? 242 ! 243 IF( ln_nicep ) CALL lim_ctl ! alerts in case of model crash 244 ! 245 ENDIF 246 #else 236 247 IF( ( MOD( kt+nn_fsbc-1, ninfo ) == 0 .OR. ntmoy == 1 ) .AND. .NOT. lk_mpp ) & 237 248 & CALL lim_dia … … 242 253 IF( ln_nicep ) CALL lim_ctl ! alerts in case of model crash 243 254 ! 255 #endif 244 256 ENDIF ! End sea-ice time step only 245 257 -
trunk/NEMO/OPA_SRC/SBC/sbcice_lim_2.F90
r1226 r1239 183 183 ! Ice model step ! 184 184 ! ---------------- ! 185 186 #if defined key_zdftke2 187 IF ( kt .NE. nitend+1 ) THEN 185 188 CALL lim_rst_opn_2 ( kt ) ! Open Ice restart file 189 ENDIF 190 #else 191 CALL lim_rst_opn_2 ( kt ) ! Open Ice restart file 192 #endif 186 193 IF( .NOT. lk_c1d ) THEN ! Ice dynamics & transport (not in 1D case) 187 194 CALL lim_dyn_2 ( kt ) ! Ice dynamics ( rheology/dynamics ) … … 197 204 CALL lim_thd_2 ( kt ) ! Ice thermodynamics 198 205 CALL lim_sbc_2 ( kt ) ! Ice/Ocean Mass & Heat fluxes 206 #if defined key_zdftke2 207 IF( kt .NE. nitend+1 ) THEN 208 IF( ( MOD( kt+nn_fsbc-1, ninfo ) == 0 .OR. ntmoy == 1 ) .AND. .NOT. lk_mpp ) & 209 & CALL lim_dia_2 ( kt ) ! Ice Diagnostics 210 CALL lim_wri_2 ( kt ) ! Ice outputs 211 IF( lrst_ice ) CALL lim_rst_write_2( kt ) ! Ice restart file 212 ! 213 ENDIF 214 #else 199 215 IF( ( MOD( kt+nn_fsbc-1, ninfo ) == 0 .OR. ntmoy == 1 ) .AND. .NOT. lk_mpp ) & 200 & CALL lim_dia_2 ( kt ) ! Ice Diagnostics 201 CALL lim_wri_2 ( kt ) ! Ice outputs 202 IF( lrst_ice ) CALL lim_rst_write_2( kt ) ! Ice restart file 216 & CALL lim_dia_2 ( kt ) ! Ice Diagnostics 217 CALL lim_wri_2 ( kt ) ! Ice outputs 218 IF( lrst_ice ) CALL lim_rst_write_2( kt ) ! Ice restart file 203 219 ! 220 #endif 204 221 ENDIF 205 222 ! -
trunk/NEMO/OPA_SRC/TRA/trazdf.F90
r1146 r1239 126 126 !!---------------------------------------------------------------------- 127 127 USE zdftke 128 USE zdftke2 128 129 USE zdfkpp 129 130 !!---------------------------------------------------------------------- … … 140 141 141 142 ! Force implicit schemes 142 IF( lk_zdftke .OR. lk_zdf kpp )nzdf = 1 ! TKE or KPP physics143 IF( ln_traldf_iso ) nzdf = 1 ! iso-neutral lateral physics144 IF( ln_traldf_hor .AND. ln_sco ) nzdf = 1 ! horizontal lateral physics in s-coordinate143 IF( lk_zdftke .OR. lk_zdftke2 .OR. lk_zdfkpp ) nzdf = 1 ! TKE or KPP physics 144 IF( ln_traldf_iso ) nzdf = 1 ! iso-neutral lateral physics 145 IF( ln_traldf_hor .AND. ln_sco ) nzdf = 1 ! horizontal lateral physics in s-coordinate 145 146 146 147 IF( ln_zdfexp .AND. nzdf == 1 ) THEN -
trunk/NEMO/OPA_SRC/ZDF/zdfini.F90
r1156 r1239 14 14 USE zdf_oce ! TKE vertical mixing 15 15 USE lib_mpp ! distribued memory computing 16 USE zdftke ! TKE vertical mixing 16 USE zdftke ! TKE vertical mixing 17 USE zdftke2 ! TKE2 vertical mixing 17 18 USE zdfkpp ! KPP vertical mixing 18 19 USE zdfddm ! double diffusion mixing … … 104 105 ioptio = ioptio+1 105 106 ENDIF 107 IF( lk_zdftke2 ) THEN 108 IF(lwp) WRITE(numout,*) ' TKE2 dependent eddy coefficients' 109 ioptio = ioptio+1 110 ENDIF 106 111 IF( lk_zdfkpp ) THEN 107 112 IF(lwp) WRITE(numout,*) ' KPP dependent eddy coefficients' … … 123 128 ioptio = ioptio+1 124 129 ENDIF 125 IF( lk_zdftke ) THEN130 IF( lk_zdftke .OR. lk_zdftke2 ) THEN 126 131 IF(lwp) WRITE(numout,*) ' use the 1.5 turbulent closure' 127 132 ENDIF … … 135 140 IF ( ioptio > 1 .AND. .NOT. lk_esopa ) & 136 141 CALL ctl_stop( ' chose between ln_zdfnpc', ' and ln_zdfevd' ) 137 IF( ioptio == 0 .AND. .NOT.( lk_zdftke .OR. lk_zdf kpp ) ) &142 IF( ioptio == 0 .AND. .NOT.( lk_zdftke .OR. lk_zdftke2 .OR. lk_zdfkpp ) ) & 138 143 CALL ctl_stop( ' except for TKE sor KPP physics, a convection scheme is', & 139 144 & ' required: ln_zdfevd or ln_zdfnpc logicals' ) -
trunk/NEMO/OPA_SRC/restart.F90
r1229 r1239 24 24 USE eosbn2 ! equation of state (eos bn2 routine) 25 25 USE trdmld_oce ! ocean active mixed layer tracers trends variables 26 #if defined key_zdftke2 27 USE zdf_oce 28 #endif 26 29 27 30 IMPLICIT NONE … … 64 67 lrst_oce = .FALSE. 65 68 nitrst = nitend 69 #if defined key_zdftke2 70 nitrst_tke2 = nitrst + 1 71 #endif 66 72 ENDIF 67 73 IF( MOD( kt - 1, nstock ) == 0 ) THEN … … 70 76 IF( nitrst > nitend ) nitrst = nitend ! make sure we write a restart at the end of the run 71 77 ENDIF 72 78 #if defined key_zdftke2 79 IF ( nitrst_tke2 .NE. kt ) nitrst_tke2 = nitrst + 1 80 #endif 73 81 ! to get better performances with NetCDF format: 74 82 ! we open and define the ocean restart file one time step before writing the data (-> at nitrst - 1) … … 110 118 INTEGER, INTENT(in) :: kt ! ocean time-step 111 119 !!---------------------------------------------------------------------- 120 121 #if defined key_zdftke2 122 IF( kt == nitrst_tke2 ) THEN 123 CALL iom_close( numrow ) ! close the restart file (only at last time step) 124 IF( .NOT. lk_trdmld ) lrst_oce = .FALSE. 125 ELSE 126 #endif 112 127 ! ! the begining of the run [s] 113 128 CALL iom_rstput( kt, nitrst, numrow, 'rdt' , rdt ) ! dynamics time step … … 129 144 CALL iom_rstput( kt, nitrst, numrow, 'hdivn' , hdivn ) 130 145 146 #if defined key_zdftke2 147 CALL iom_rstput( kt, nitrst, numrow, 'rhop' , rhop ) 148 #endif 131 149 IF( nn_dynhpg_rst == 1 .OR. lk_vvl ) THEN 132 150 CALL iom_rstput( kt, nitrst, numrow, 'rhd' , rhd ) 151 #if defined key_zdftke2 152 CALL iom_rstput( kt, nitrst, numrow, 'rn2' , rn2 ) 153 CALL iom_rstput( kt, nitrst, numrow, 'avt' , avt ) 154 ENDIF 155 #else 133 156 CALL iom_rstput( kt, nitrst, numrow, 'rhop', rhop ) 157 #endif 134 158 IF( ln_zps ) THEN 135 159 CALL iom_rstput( kt, nitrst, numrow, 'gtu' , gtu ) … … 141 165 ENDIF 142 166 ENDIF 143 167 #if ! defined key_zdftke2 144 168 IF( kt == nitrst ) THEN 145 169 CALL iom_close( numrow ) ! close the restart file (only at last time step) 146 170 IF( .NOT. lk_trdmld ) lrst_oce = .FALSE. 147 171 ENDIF 172 #endif 148 173 ! 149 174 END SUBROUTINE rst_write … … 163 188 !! - barotropic stream function arrays ("key_dynspg_rl" defined) 164 189 !! or free surface arrays 165 !! - tke arrays (lk_zdftke=T )190 !! - tke arrays (lk_zdftke=T .OR. lk_zdftke2=T) 166 191 !! for this last three records, the previous characteristics 167 192 !! could be different with those used in the present run. … … 215 240 ENDIF 216 241 242 #if defined key_zdftke2 243 CALL eos( tb, sb, rhd, rhop ) ! before potential and in situ densities 244 IF( iom_varid( numror, 'rhd', ldstop = .FALSE. ) > 0 ) THEN 245 CALL iom_get( numror, jpdom_autoglo, 'rhd' , rhd ) 246 ENDIF 247 IF( iom_varid( numror, 'rhop', ldstop = .FALSE. ) > 0 ) THEN 248 CALL iom_get( numror, jpdom_autoglo, 'rhop', rhop ) 249 ENDIF 250 #else 217 251 IF( iom_varid( numror, 'rhd', ldstop = .FALSE. ) > 0 ) THEN 218 252 CALL iom_get( numror, jpdom_autoglo, 'rhd' , rhd ) … … 221 255 CALL eos( tb, sb, rhd, rhop ) ! before potential and in situ densities 222 256 ENDIF 257 #endif 258 #if defined key_zdftke2 259 CALL eos_init ! usefull to get the equation state type neos parameter 260 IF( iom_varid( numror, 'rn2', ldstop = .FALSE. ) > 0 ) THEN 261 CALL iom_get( numror, jpdom_autoglo, 'rn2' , rn2 ) 262 ELSE 263 IF ( ln_dynhpg_imp ) THEN 264 CALL bn2( tb, sb, rn2 ) ! before Brunt-Vaisala frequency 265 ENDIF 266 ENDIF 267 IF( iom_varid( numror, 'avt', ldstop = .FALSE. ) > 0 ) THEN 268 CALL iom_get( numror, jpdom_autoglo, 'avt' , avt ) 269 ELSE 270 IF ( ln_dynhpg_imp ) avt (:,:,:) = 1.2e-5 * tmask(:,:,:) 271 ENDIF 272 #endif 273 223 274 IF( ln_zps .AND. .NOT. lk_c1d ) THEN 224 275 IF( iom_varid( numror, 'gtu', ldstop = .FALSE. ) > 0 ) THEN -
trunk/NEMO/OPA_SRC/step.F90
r1218 r1239 88 88 USE zdfbfr ! bottom friction (zdf_bfr routine) 89 89 USE zdftke ! TKE vertical mixing (zdf_tke routine) 90 USE zdftke2 ! TKE2 vertical mixing (zdf_tke2 routine) 90 91 USE zdfkpp ! KPP vertical mixing (zdf_kpp routine) 91 92 USE zdfddm ! double diffusion mixing (zdf_ddm routine) … … 202 203 ! Ocean physics update 203 204 !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 205 #if defined key_zdftke2 206 IF ( ln_dynhpg_imp ) THEN 207 !----------------------------------------------------------------------- 208 ! LATERAL PHYSICS 209 !----------------------------------------------------------------------- 210 ! N.B. ua, va, ta, sa arrays are used as workspace in this section 211 !----------------------------------------------------------------------- 212 CALL zdf_mxl( kstp ) ! mixed layer depth 213 IF( lk_ldfslp ) CALL ldf_slp( kstp, rhd, rn2 ) ! before slope of the lateral mixing 214 # if defined key_traldf_c2d 215 IF( lk_traldf_eiv ) CALL ldf_eiv( kstp ) ! eddy induced velocity coefficient 216 # endif 217 ENDIF 218 #endif 204 219 !----------------------------------------------------------------------- 205 220 ! VERTICAL PHYSICS … … 207 222 ! N.B. ua, va, ta, sa arrays are used as workspace in this section 208 223 !----------------------------------------------------------------------- 209 224 #if defined key_zdftke2 225 CALL bn2( tn, sn, rn2 ) ! now Brunt-Vaisala frequency 226 #else 210 227 CALL bn2( tb, sb, rn2 ) ! before Brunt-Vaisala frequency 211 228 #endif 212 229 ! ! Vertical eddy viscosity and diffusivity coefficients 213 230 IF( lk_zdfric ) CALL zdf_ric( kstp ) ! Richardson number dependent Kz 214 231 215 IF( lk_zdftke ) CALL zdf_tke( kstp ) ! TKE closure scheme for Kz 232 IF( lk_zdftke ) CALL zdf_tke ( kstp ) ! TKE closure scheme for Kz 233 IF( lk_zdftke2) CALL zdf_tke2( kstp ) ! TKE2 closure scheme for Kz 216 234 217 235 IF( lk_zdfkpp ) CALL zdf_kpp( kstp ) ! KPP closure scheme for Kz … … 236 254 CALL zdf_mxl( kstp ) ! mixed layer depth 237 255 238 256 #if defined key_zdftke2 257 IF( .NOT. ln_dynhpg_imp ) THEN 258 CALL bn2( tb, sb, rn2 ) ! before Brunt-Vaisala frequency 259 CALL eos( tb, sb, rhd, rhop ) ! now (swap=before) in situ density for dynhpg module 260 #endif 239 261 !----------------------------------------------------------------------- 240 262 ! LATERAL PHYSICS … … 245 267 #if defined key_traldf_c2d 246 268 IF( lk_traldf_eiv ) CALL ldf_eiv( kstp ) ! eddy induced velocity coefficient 269 # endif 270 #if defined key_zdftke2 271 ENDIF 247 272 #endif 248 273 … … 278 303 #endif 279 304 CALL tra_zdf ( kstp ) ! vertical mixing and after tracer fields 305 306 #if ! defined key_zdftke2 280 307 IF( ln_zdfnpc ) CALL tra_npc ( kstp ) ! update after fields by non-penetrative convection 281 308 CALL tra_nxt ( kstp ) ! tracer fields at next time step … … 292 319 & gtv, gsv, grv ) 293 320 ENDIF 294 321 #else 322 IF( .NOT. ln_dynhpg_imp ) THEN ! centered hpg (default case) 323 CALL eos( tn, sn, rhd, rhop ) ! now (swap=before) in situ density for dynhpg module 324 IF( ln_zps ) CALL zps_hde( kstp, tn, sn, rhd, & ! Partial steps: now horizontal gradient 325 & gtu, gsu, gru, & ! of t, s, rd at the bottom ocean level 326 & gtv, gsv, grv ) 327 ENDIF 328 IF( ln_zdfnpc ) CALL tra_npc ( kstp ) ! update after fields by non-penetrative convection 329 CALL tra_nxt ( kstp ) ! tracer fields at next time step 330 IF( ln_dynhpg_imp ) THEN ! semi-implicit hpg 331 CALL eos( ta, sa, rhd, rhop ) ! Time-filtered in situ density used in dynhpg module 332 IF( lk_ldfslp ) CALL bn2( ta, sa, rn2 ) ! Time-filtered Brunt-Vaisala frequency 333 IF( ln_zps ) CALL zps_hde( kstp, ta, sa, rhd, & ! Partial steps: time filtered hor. gradient 334 & gtu, gsu, gru, & ! of t, s, rd at the bottom ocean level 335 & gtv, gsv, grv ) 336 ENDIF 337 #endif 295 338 !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 296 339 ! Dynamics … … 363 406 364 407 IF( lk_cpl ) CALL sbc_cpl_snd( kstp ) ! coupled mode : field exchanges 408 409 #if defined key_zdftke2 410 IF( ( kstp == nitend ).AND.( lrst_oce ) ) THEN 411 412 CALL day( kstp+1 ) ! Calendar 413 !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 414 ! Update data, open boundaries, surface boundary condition (including sea-ice) 415 !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 416 CALL sbc ( kstp+1 ) ! Sea Boundary Condition (including sea-ice) 417 !----------------------------------------------------------------------- 418 ! VERTICAL PHYSICS 419 !----------------------------------------------------------------------- 420 ! N.B. ua, va, ta, sa arrays are used as workspace in this section 421 !----------------------------------------------------------------------- 422 CALL bn2( tn, sn, rn2 ) ! now Brunt-Vaisala frequency 423 ! ! Vertical eddy viscosity and diffusivity coefficients 424 IF( lk_zdftke2 ) CALL zdf_tke2 ( kstp+1 ) ! TKE2 closure scheme for Kz 425 CALL rst_write( kstp+1 ) ! close the restart file 426 ENDIF 427 #endif 365 428 ! 366 429 END SUBROUTINE stp -
trunk/NEMO/TOP_SRC/TRP/trctrp_ctl.F90
r1146 r1239 249 249 ! ... vertical mixing 250 250 ! time stepping scheme (N.B. TKE && KPP scheme => force the use of implicit scheme) 251 #if defined key_zdftke || defined key_zdf kpp251 #if defined key_zdftke || defined key_zdftke2 || defined key_zdfkpp 252 252 l_trczdf_exp = .FALSE. ! use implicit scheme 253 253 l_trczdf_imp = .TRUE.
Note: See TracChangeset
for help on using the changeset viewer.