[3] | 1 | MODULE istate |
---|
| 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE istate *** |
---|
| 4 | !! Ocean state : initial state setting |
---|
| 5 | !!===================================================================== |
---|
[508] | 6 | !! History : 4.0 ! 89-12 (P. Andrich) Original code |
---|
| 7 | !! 5.0 ! 91-11 (G. Madec) rewritting |
---|
| 8 | !! 6.0 ! 96-01 (G. Madec) terrain following coordinates |
---|
| 9 | !! 8.0 ! 01-09 (M. Levy, M. Ben Jelloul) istate_eel |
---|
| 10 | !! 8.0 ! 01-09 (M. Levy, M. Ben Jelloul) istate_uvg |
---|
| 11 | !! 9.0 ! 03-08 (G. Madec) F90: Free form, modules |
---|
| 12 | !! 9.0 ! 03-09 (G. Madec, C. Talandier) add EEL R5 |
---|
| 13 | !! 9.0 ! 04-05 (A. Koch-Larrouy) istate_gyre |
---|
| 14 | !! 9.0 ! 06-07 (S. Masson) distributed restart using iom |
---|
| 15 | !!---------------------------------------------------------------------- |
---|
[3] | 16 | |
---|
| 17 | !!---------------------------------------------------------------------- |
---|
| 18 | !! istate_init : initial state setting |
---|
| 19 | !! istate_tem : analytical profile for initial Temperature |
---|
| 20 | !! istate_sal : analytical profile for initial Salinity |
---|
| 21 | !! istate_eel : initial state setting of EEL R5 configuration |
---|
[93] | 22 | !! istate_gyre : initial state setting of GYRE configuration |
---|
[3] | 23 | !! istate_uvg : initial velocity in geostropic balance |
---|
| 24 | !!---------------------------------------------------------------------- |
---|
| 25 | USE oce ! ocean dynamics and active tracers |
---|
| 26 | USE dom_oce ! ocean space and time domain |
---|
| 27 | USE daymod ! |
---|
| 28 | USE ldftra_oce ! ocean active tracers: lateral physics |
---|
| 29 | USE zdf_oce ! ocean vertical physics |
---|
| 30 | USE phycst ! physical constants |
---|
| 31 | USE dtatem ! temperature data (dta_tem routine) |
---|
| 32 | USE dtasal ! salinity data (dta_sal routine) |
---|
| 33 | USE restart ! ocean restart (rst_read routine) |
---|
[508] | 34 | USE in_out_manager ! I/O manager |
---|
| 35 | USE iom |
---|
[900] | 36 | USE c1d ! re-initialization of u-v mask for the 1D configuration |
---|
[544] | 37 | USE zpshde ! partial step: hor. derivative (zps_hde routine) |
---|
| 38 | USE eosbn2 ! equation of state (eos bn2 routine) |
---|
[593] | 39 | USE domvvl ! varying vertical mesh |
---|
| 40 | USE dynspg_oce ! pressure gradient schemes |
---|
| 41 | USE dynspg_flt ! pressure gradient schemes |
---|
| 42 | USE dynspg_exp ! pressure gradient schemes |
---|
| 43 | USE dynspg_ts ! pressure gradient schemes |
---|
[508] | 44 | |
---|
[3] | 45 | IMPLICIT NONE |
---|
| 46 | PRIVATE |
---|
| 47 | |
---|
[508] | 48 | PUBLIC istate_init ! routine called by step.F90 |
---|
[3] | 49 | |
---|
| 50 | !! * Substitutions |
---|
| 51 | # include "domzgr_substitute.h90" |
---|
| 52 | # include "vectopt_loop_substitute.h90" |
---|
| 53 | !!---------------------------------------------------------------------- |
---|
[508] | 54 | !! OPA 9.0 , LOCEAN-IPSL (2006) |
---|
[888] | 55 | !! $Id$ |
---|
[508] | 56 | !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt) |
---|
[3] | 57 | !!---------------------------------------------------------------------- |
---|
| 58 | |
---|
| 59 | CONTAINS |
---|
| 60 | |
---|
| 61 | SUBROUTINE istate_init |
---|
| 62 | !!---------------------------------------------------------------------- |
---|
| 63 | !! *** ROUTINE istate_init *** |
---|
| 64 | !! |
---|
[508] | 65 | !! ** Purpose : Initialization of the dynamics and tracer fields. |
---|
[3] | 66 | !!---------------------------------------------------------------------- |
---|
[544] | 67 | USE eosbn2 ! eq. of state, Brunt Vaisala frequency (eos routine) |
---|
[2068] | 68 | ! - ML - needed for initialization of e3t_b |
---|
| 69 | INTEGER :: jk ! dummy loop indice |
---|
[3] | 70 | |
---|
[508] | 71 | IF(lwp) WRITE(numout,*) |
---|
| 72 | IF(lwp) WRITE(numout,*) 'istate_ini : Initialization of the dynamics and tracers' |
---|
| 73 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~' |
---|
[3] | 74 | |
---|
| 75 | rhd (:,:,:) = 0.e0 |
---|
| 76 | rhop (:,:,:) = 0.e0 |
---|
| 77 | rn2 (:,:,:) = 0.e0 |
---|
| 78 | |
---|
[15] | 79 | IF( ln_rstart ) THEN ! Restart from a file |
---|
[3] | 80 | ! ! ------------------- |
---|
| 81 | neuler = 1 ! Set time-step indicator at nit000 (leap-frog) |
---|
| 82 | CALL rst_read ! Read the restart file |
---|
[1130] | 83 | CALL day_init ! model calendar (using both namelist and restart infos) |
---|
[3] | 84 | ELSE |
---|
| 85 | ! ! Start from rest |
---|
| 86 | ! ! --------------- |
---|
[1130] | 87 | numror = 0 ! define numror = 0 -> no restart file to read |
---|
[3] | 88 | neuler = 0 ! Set time-step indicator at nit000 (euler forward) |
---|
[1130] | 89 | CALL day_init ! model calendar (using both namelist and restart infos) |
---|
[508] | 90 | ! ! Initialization of ocean to zero |
---|
| 91 | ! before fields ! now fields |
---|
[684] | 92 | ub (:,:,:) = 0.e0 ; un (:,:,:) = 0.e0 ; sshb(:,:) = 0.e0 |
---|
| 93 | vb (:,:,:) = 0.e0 ; vn (:,:,:) = 0.e0 ; sshn(:,:) = 0.e0 |
---|
| 94 | rotb (:,:,:) = 0.e0 ; rotn (:,:,:) = 0.e0 |
---|
| 95 | hdivb(:,:,:) = 0.e0 ; hdivn(:,:,:) = 0.e0 |
---|
[508] | 96 | ! |
---|
[3] | 97 | IF( cp_cfg == 'eel' ) THEN |
---|
[93] | 98 | CALL istate_eel ! EEL configuration : start from pre-defined |
---|
| 99 | ! ! velocity and thermohaline fields |
---|
[434] | 100 | ELSEIF( cp_cfg == 'gyre' ) THEN |
---|
[93] | 101 | CALL istate_gyre ! GYRE configuration : start from pre-defined temperature |
---|
| 102 | ! ! and salinity fields |
---|
[3] | 103 | ELSE |
---|
[508] | 104 | ! ! Other configurations: Initial temperature and salinity fields |
---|
[3] | 105 | #if defined key_dtatem |
---|
| 106 | CALL dta_tem( nit000 ) ! read 3D temperature data |
---|
| 107 | tb(:,:,:) = t_dta(:,:,:) ! use temperature data read |
---|
| 108 | tn(:,:,:) = t_dta(:,:,:) |
---|
| 109 | #else |
---|
| 110 | IF(lwp) WRITE(numout,*) ! analytical temperature profile |
---|
[508] | 111 | IF(lwp) WRITE(numout,*)' Temperature initialization using an analytic profile' |
---|
[3] | 112 | CALL istate_tem |
---|
| 113 | #endif |
---|
| 114 | #if defined key_dtasal |
---|
| 115 | CALL dta_sal( nit000 ) ! read 3D salinity data |
---|
| 116 | sb(:,:,:) = s_dta(:,:,:) ! use salinity data read |
---|
| 117 | sn(:,:,:) = s_dta(:,:,:) |
---|
| 118 | #else |
---|
| 119 | ! No salinity data |
---|
| 120 | IF(lwp)WRITE(numout,*) ! analytical salinity profile |
---|
[508] | 121 | IF(lwp)WRITE(numout,*)' Salinity initialisation using a constant value' |
---|
[3] | 122 | CALL istate_sal |
---|
| 123 | #endif |
---|
| 124 | ENDIF |
---|
| 125 | |
---|
[544] | 126 | CALL eos( tb, sb, rhd, rhop ) ! before potential and in situ densities |
---|
| 127 | |
---|
[899] | 128 | IF( ln_zps .AND. .NOT. lk_c1d ) & |
---|
[544] | 129 | & CALL zps_hde( nit000, tb, sb, rhd, & ! Partial steps: before Horizontal DErivative |
---|
| 130 | & gtu, gsu, gru, & ! of t, s, rd at the bottom ocean level |
---|
| 131 | & gtv, gsv, grv ) |
---|
[2068] | 132 | |
---|
| 133 | ! - ML - sshn could be modified by istate_eel, so that initialization of fse3t_b is done here |
---|
| 134 | IF( lk_vvl ) THEN |
---|
| 135 | DO jk = 1, jpk |
---|
| 136 | fse3t_b(:,:,jk) = fse3t_n(:,:,jk) |
---|
| 137 | ENDDO |
---|
| 138 | ENDIF |
---|
[3] | 139 | ENDIF |
---|
[1438] | 140 | ! |
---|
| 141 | IF( lk_agrif ) THEN |
---|
[593] | 142 | ! read free surface arrays in restart file |
---|
| 143 | IF( ln_rstart ) THEN |
---|
| 144 | IF( lk_dynspg_flt ) CALL flt_rst( nit000, 'READ' ) ! read or initialize the following fields |
---|
[1438] | 145 | ! ! gcx, gcxb for agrif_opa_init |
---|
| 146 | ENDIF ! explicit case not coded yet with AGRIF |
---|
[1200] | 147 | ENDIF |
---|
[508] | 148 | ! |
---|
[3] | 149 | END SUBROUTINE istate_init |
---|
| 150 | |
---|
| 151 | |
---|
| 152 | SUBROUTINE istate_tem |
---|
| 153 | !!--------------------------------------------------------------------- |
---|
| 154 | !! *** ROUTINE istate_tem *** |
---|
| 155 | !! |
---|
| 156 | !! ** Purpose : Intialization of the temperature field with an |
---|
| 157 | !! analytical profile or a file (i.e. in EEL configuration) |
---|
| 158 | !! |
---|
| 159 | !! ** Method : Use Philander analytic profile of temperature |
---|
| 160 | !! |
---|
| 161 | !! References : Philander ??? |
---|
| 162 | !!---------------------------------------------------------------------- |
---|
| 163 | INTEGER :: ji, jj, jk |
---|
| 164 | !!---------------------------------------------------------------------- |
---|
[508] | 165 | ! |
---|
[3] | 166 | IF(lwp) WRITE(numout,*) |
---|
| 167 | IF(lwp) WRITE(numout,*) 'istate_tem : initial temperature profile' |
---|
| 168 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~' |
---|
| 169 | |
---|
| 170 | DO jk = 1, jpk |
---|
| 171 | DO jj = 1, jpj |
---|
| 172 | DO ji = 1, jpi |
---|
| 173 | tn(ji,jj,jk) = ( ( ( 7.5 - 0.*ABS(gphit(ji,jj))/30. ) & |
---|
[15] | 174 | & *( 1.-TANH((fsdept(ji,jj,jk)-80.)/30.) ) & |
---|
| 175 | & + 10.*(5000.-fsdept(ji,jj,jk))/5000.) ) * tmask(ji,jj,jk) |
---|
[3] | 176 | tb(ji,jj,jk) = tn(ji,jj,jk) |
---|
| 177 | END DO |
---|
| 178 | END DO |
---|
| 179 | END DO |
---|
| 180 | |
---|
[79] | 181 | IF(lwp) CALL prizre( tn , jpi , jpj , jpk , jpj/2 , & |
---|
| 182 | & 1 , jpi , 5 , 1 , jpk , & |
---|
| 183 | & 1 , 1. , numout ) |
---|
[508] | 184 | ! |
---|
[3] | 185 | END SUBROUTINE istate_tem |
---|
| 186 | |
---|
| 187 | |
---|
| 188 | SUBROUTINE istate_sal |
---|
| 189 | !!--------------------------------------------------------------------- |
---|
| 190 | !! *** ROUTINE istate_sal *** |
---|
| 191 | !! |
---|
| 192 | !! ** Purpose : Intialize the salinity field with an analytic profile |
---|
| 193 | !! |
---|
| 194 | !! ** Method : Use to a constant value 35.5 |
---|
| 195 | !! |
---|
| 196 | !! ** Action : Initialize sn and sb |
---|
| 197 | !!---------------------------------------------------------------------- |
---|
| 198 | REAL(wp) :: zsal = 35.50_wp |
---|
| 199 | !!---------------------------------------------------------------------- |
---|
| 200 | |
---|
| 201 | IF(lwp) WRITE(numout,*) |
---|
| 202 | IF(lwp) WRITE(numout,*) 'istate_sal : initial salinity : ', zsal |
---|
| 203 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~' |
---|
| 204 | |
---|
| 205 | sn(:,:,:) = zsal * tmask(:,:,:) |
---|
| 206 | sb(:,:,:) = sn(:,:,:) |
---|
| 207 | |
---|
| 208 | END SUBROUTINE istate_sal |
---|
| 209 | |
---|
| 210 | |
---|
| 211 | SUBROUTINE istate_eel |
---|
| 212 | !!---------------------------------------------------------------------- |
---|
| 213 | !! *** ROUTINE istate_eel *** |
---|
| 214 | !! |
---|
| 215 | !! ** Purpose : Initialization of the dynamics and tracers for EEL R5 |
---|
| 216 | !! configuration (channel with or without a topographic bump) |
---|
| 217 | !! |
---|
| 218 | !! ** Method : - set temprature field |
---|
| 219 | !! - set salinity field |
---|
| 220 | !! - set velocity field including horizontal divergence |
---|
| 221 | !! and relative vorticity fields |
---|
| 222 | !!---------------------------------------------------------------------- |
---|
| 223 | USE eosbn2 ! eq. of state, Brunt Vaisala frequency (eos routine) |
---|
| 224 | USE divcur ! hor. divergence & rel. vorticity (div_cur routine) |
---|
[473] | 225 | USE iom |
---|
[3] | 226 | |
---|
| 227 | INTEGER :: inum ! temporary logical unit |
---|
| 228 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
[479] | 229 | INTEGER :: ijloc |
---|
[508] | 230 | REAL(wp) :: zh1, zh2, zslope, zcst, zfcor ! temporary scalars |
---|
[593] | 231 | REAL(wp) :: zt1 = 15._wp, & ! surface temperature value (EEL R5) |
---|
| 232 | & zt2 = 5._wp, & ! bottom temperature value (EEL R5) |
---|
| 233 | & zsal = 35.0_wp, & ! constant salinity (EEL R2, R5 and R6) |
---|
[508] | 234 | & zueel = 0.1_wp ! constant uniform zonal velocity (EEL R5) |
---|
| 235 | REAL(wp), DIMENSION(jpiglo,jpjglo) :: zssh ! initial ssh over the global domain |
---|
[3] | 236 | !!---------------------------------------------------------------------- |
---|
| 237 | |
---|
| 238 | SELECT CASE ( jp_cfg ) |
---|
| 239 | ! ! ==================== |
---|
| 240 | CASE ( 5 ) ! EEL R5 configuration |
---|
| 241 | ! ! ==================== |
---|
| 242 | |
---|
| 243 | ! set temperature field with a linear profile |
---|
| 244 | ! ------------------------------------------- |
---|
| 245 | IF(lwp) WRITE(numout,*) |
---|
| 246 | IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: linear temperature profile' |
---|
| 247 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~' |
---|
| 248 | |
---|
[467] | 249 | zh1 = gdept_0( 1 ) |
---|
| 250 | zh2 = gdept_0(jpkm1) |
---|
[3] | 251 | |
---|
| 252 | zslope = ( zt1 - zt2 ) / ( zh1 - zh2 ) |
---|
| 253 | zcst = ( zt1 * ( zh1 - zh2) - ( zt1 - zt2 ) * zh1 ) / ( zh1 - zh2 ) |
---|
| 254 | |
---|
| 255 | DO jk = 1, jpk |
---|
[593] | 256 | tn(:,:,jk) = ( zt2 + zt1 * exp( - fsdept(:,:,jk) / 1000 ) ) * tmask(:,:,jk) |
---|
[3] | 257 | tb(:,:,jk) = tn(:,:,jk) |
---|
| 258 | END DO |
---|
| 259 | |
---|
| 260 | IF(lwp) CALL prizre( tn , jpi , jpj , jpk , jpj/2 , & |
---|
| 261 | & 1 , jpi , 5 , 1 , jpk , & |
---|
| 262 | & 1 , 1. , numout ) |
---|
| 263 | |
---|
| 264 | ! set salinity field to a constant value |
---|
| 265 | ! -------------------------------------- |
---|
| 266 | IF(lwp) WRITE(numout,*) |
---|
| 267 | IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal |
---|
| 268 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~' |
---|
| 269 | |
---|
| 270 | sn(:,:,:) = zsal * tmask(:,:,:) |
---|
| 271 | sb(:,:,:) = sn(:,:,:) |
---|
| 272 | |
---|
| 273 | |
---|
| 274 | ! set the dynamics: U,V, hdiv, rot (and ssh if necessary) |
---|
| 275 | ! ---------------- |
---|
| 276 | ! Start EEL5 configuration with barotropic geostrophic velocities |
---|
| 277 | ! according the sshb and sshn SSH imposed. |
---|
[479] | 278 | ! we assume a uniform grid (hence the use of e1t(1,1) for delta_y) |
---|
| 279 | ! we use the Coriolis frequency at mid-channel. |
---|
| 280 | |
---|
| 281 | ub(:,:,:) = zueel * umask(:,:,:) |
---|
[3] | 282 | un(:,:,:) = ub(:,:,:) |
---|
[479] | 283 | ijloc = mj0(INT(jpjglo-1)/2) |
---|
| 284 | zfcor = ff(1,ijloc) |
---|
[3] | 285 | |
---|
| 286 | DO jj = 1, jpjglo |
---|
[479] | 287 | zssh(:,jj) = - (FLOAT(jj)- FLOAT(jpjglo-1)/2.)*zueel*e1t(1,1)*zfcor/grav |
---|
[3] | 288 | END DO |
---|
[479] | 289 | |
---|
| 290 | IF(lwp) THEN |
---|
| 291 | WRITE(numout,*) ' Uniform zonal velocity for EEL R5:',zueel |
---|
| 292 | WRITE(numout,*) ' Geostrophic SSH profile as a function of y:' |
---|
| 293 | WRITE(numout,'(12(1x,f6.2))') zssh(1,:) |
---|
| 294 | ENDIF |
---|
| 295 | |
---|
[3] | 296 | DO jj = 1, nlcj |
---|
| 297 | DO ji = 1, nlci |
---|
| 298 | sshb(ji,jj) = zssh( mig(ji) , mjg(jj) ) * tmask(ji,jj,1) |
---|
| 299 | END DO |
---|
| 300 | END DO |
---|
| 301 | sshb(nlci+1:jpi, : ) = 0.e0 ! set to zero extra mpp columns |
---|
| 302 | sshb( : ,nlcj+1:jpj) = 0.e0 ! set to zero extra mpp rows |
---|
| 303 | |
---|
| 304 | sshn(:,:) = sshb(:,:) ! set now ssh to the before value |
---|
| 305 | |
---|
[593] | 306 | IF( nn_rstssh /= 0 ) THEN |
---|
| 307 | nn_rstssh = 0 ! hand-made initilization of ssh |
---|
| 308 | CALL ctl_warn( 'istate_eel: force nn_rstssh = 0' ) |
---|
[558] | 309 | ENDIF |
---|
| 310 | |
---|
[3] | 311 | ! horizontal divergence and relative vorticity (curl) |
---|
| 312 | CALL div_cur( nit000 ) |
---|
| 313 | |
---|
| 314 | ! N.B. the vertical velocity will be computed from the horizontal divergence field |
---|
| 315 | ! in istate by a call to wzv routine |
---|
| 316 | |
---|
| 317 | |
---|
| 318 | ! ! ========================== |
---|
| 319 | CASE ( 2 , 6 ) ! EEL R2 or R6 configuration |
---|
| 320 | ! ! ========================== |
---|
| 321 | |
---|
| 322 | ! set temperature field with a NetCDF file |
---|
| 323 | ! ---------------------------------------- |
---|
| 324 | IF(lwp) WRITE(numout,*) |
---|
| 325 | IF(lwp) WRITE(numout,*) 'istate_eel : EEL R2 or R6: read initial temperature in a NetCDF file' |
---|
| 326 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~' |
---|
| 327 | |
---|
[473] | 328 | CALL iom_open ( 'eel.initemp', inum ) |
---|
| 329 | CALL iom_get ( inum, jpdom_data, 'initemp', tb ) ! read before temprature (tb) |
---|
| 330 | CALL iom_close( inum ) |
---|
| 331 | |
---|
| 332 | tn(:,:,:) = tb(:,:,:) ! set nox temperature to tb |
---|
[3] | 333 | |
---|
| 334 | IF(lwp) CALL prizre( tn , jpi , jpj , jpk , jpj/2 , & |
---|
| 335 | & 1 , jpi , 5 , 1 , jpk , & |
---|
| 336 | & 1 , 1. , numout ) |
---|
| 337 | |
---|
| 338 | |
---|
| 339 | ! set salinity field to a constant value |
---|
| 340 | ! -------------------------------------- |
---|
| 341 | IF(lwp) WRITE(numout,*) |
---|
| 342 | IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal |
---|
| 343 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~' |
---|
| 344 | |
---|
| 345 | sn(:,:,:) = zsal * tmask(:,:,:) |
---|
| 346 | sb(:,:,:) = sn(:,:,:) |
---|
| 347 | |
---|
| 348 | ! ! =========================== |
---|
| 349 | CASE DEFAULT ! NONE existing configuration |
---|
| 350 | ! ! =========================== |
---|
[473] | 351 | WRITE(ctmp1,*) 'EEL with a ', jp_cfg,' km resolution is not coded' |
---|
| 352 | CALL ctl_stop( ctmp1 ) |
---|
| 353 | |
---|
[3] | 354 | END SELECT |
---|
| 355 | |
---|
| 356 | END SUBROUTINE istate_eel |
---|
| 357 | |
---|
| 358 | |
---|
[93] | 359 | SUBROUTINE istate_gyre |
---|
| 360 | !!---------------------------------------------------------------------- |
---|
| 361 | !! *** ROUTINE istate_gyre *** |
---|
| 362 | !! |
---|
| 363 | !! ** Purpose : Initialization of the dynamics and tracers for GYRE |
---|
| 364 | !! configuration (double gyre with rotated domain) |
---|
| 365 | !! |
---|
| 366 | !! ** Method : - set temprature field |
---|
| 367 | !! - set salinity field |
---|
| 368 | !!---------------------------------------------------------------------- |
---|
[473] | 369 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
[508] | 370 | INTEGER :: inum ! temporary logical unit |
---|
| 371 | INTEGER, PARAMETER :: ntsinit = 0 ! (0/1) (analytical/input data files) T&S initialization |
---|
[93] | 372 | !!---------------------------------------------------------------------- |
---|
| 373 | |
---|
[434] | 374 | SELECT CASE ( ntsinit) |
---|
[93] | 375 | |
---|
[434] | 376 | CASE ( 0 ) ! analytical T/S profil deduced from LEVITUS |
---|
| 377 | IF(lwp) WRITE(numout,*) |
---|
| 378 | IF(lwp) WRITE(numout,*) 'istate_gyre : initial analytical T and S profil deduced from LEVITUS ' |
---|
| 379 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~~' |
---|
[93] | 380 | |
---|
[434] | 381 | DO jk = 1, jpk |
---|
| 382 | DO jj = 1, jpj |
---|
| 383 | DO ji = 1, jpi |
---|
| 384 | tn(ji,jj,jk) = ( 16. - 12. * TANH( (fsdept(ji,jj,jk) - 400) / 700 ) ) & |
---|
| 385 | & * (-TANH( (500-fsdept(ji,jj,jk)) / 150 ) + 1) / 2 & |
---|
| 386 | & + ( 15. * ( 1. - TANH( (fsdept(ji,jj,jk)-50.) / 1500.) ) & |
---|
| 387 | & - 1.4 * TANH((fsdept(ji,jj,jk)-100.) / 100.) & |
---|
| 388 | & + 7. * (1500. - fsdept(ji,jj,jk)) / 1500. ) & |
---|
| 389 | & * (-TANH( (fsdept(ji,jj,jk) - 500) / 150) + 1) / 2 |
---|
| 390 | tn(ji,jj,jk) = tn(ji,jj,jk) * tmask(ji,jj,jk) |
---|
| 391 | tb(ji,jj,jk) = tn(ji,jj,jk) |
---|
| 392 | |
---|
| 393 | sn(ji,jj,jk) = ( 36.25 - 1.13 * TANH( (fsdept(ji,jj,jk) - 305) / 460 ) ) & |
---|
| 394 | & * (-TANH((500 - fsdept(ji,jj,jk)) / 150) + 1) / 2 & |
---|
| 395 | & + ( 35.55 + 1.25 * (5000. - fsdept(ji,jj,jk)) / 5000. & |
---|
| 396 | & - 1.62 * TANH( (fsdept(ji,jj,jk) - 60. ) / 650. ) & |
---|
| 397 | & + 0.2 * TANH( (fsdept(ji,jj,jk) - 35. ) / 100. ) & |
---|
| 398 | & + 0.2 * TANH( (fsdept(ji,jj,jk) - 1000.) / 5000.) ) & |
---|
| 399 | & * (-TANH((fsdept(ji,jj,jk) - 500) / 150) + 1) / 2 |
---|
| 400 | sn(ji,jj,jk) = sn(ji,jj,jk) * tmask(ji,jj,jk) |
---|
| 401 | sb(ji,jj,jk) = sn(ji,jj,jk) |
---|
| 402 | END DO |
---|
[93] | 403 | END DO |
---|
| 404 | END DO |
---|
| 405 | |
---|
[434] | 406 | CASE ( 1 ) ! T/S data fields read in dta_tem.nc/data_sal.nc files |
---|
| 407 | IF(lwp) WRITE(numout,*) |
---|
| 408 | IF(lwp) WRITE(numout,*) 'istate_gyre : initial T and S read from dta_tem.nc/data_sal.nc files' |
---|
| 409 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~~' |
---|
| 410 | IF(lwp) WRITE(numout,*) ' NetCDF FORMAT' |
---|
| 411 | |
---|
| 412 | ! Read temperature field |
---|
| 413 | ! ---------------------- |
---|
[473] | 414 | CALL iom_open ( 'data_tem', inum ) |
---|
| 415 | CALL iom_get ( inum, jpdom_data, 'votemper', tn ) |
---|
| 416 | CALL iom_close( inum ) |
---|
[434] | 417 | |
---|
[473] | 418 | tn(:,:,:) = tn(:,:,:) * tmask(:,:,:) |
---|
| 419 | tb(:,:,:) = tn(:,:,:) |
---|
[434] | 420 | |
---|
| 421 | ! Read salinity field |
---|
| 422 | ! ------------------- |
---|
[473] | 423 | CALL iom_open ( 'data_sal', inum ) |
---|
| 424 | CALL iom_get ( inum, jpdom_data, 'vosaline', sn ) |
---|
| 425 | CALL iom_close( inum ) |
---|
[434] | 426 | |
---|
[473] | 427 | sn(:,:,:) = sn(:,:,:) * tmask(:,:,:) |
---|
| 428 | sb(:,:,:) = sn(:,:,:) |
---|
[434] | 429 | |
---|
| 430 | END SELECT |
---|
| 431 | |
---|
[93] | 432 | IF(lwp) THEN |
---|
| 433 | WRITE(numout,*) |
---|
| 434 | WRITE(numout,*) ' Initial temperature and salinity profiles:' |
---|
[467] | 435 | WRITE(numout, "(9x,' level gdept_0 temperature salinity ')" ) |
---|
| 436 | WRITE(numout, "(10x, i4, 3f10.2)" ) ( jk, gdept_0(jk), tn(2,2,jk), sn(2,2,jk), jk = 1, jpk ) |
---|
[93] | 437 | ENDIF |
---|
| 438 | |
---|
| 439 | END SUBROUTINE istate_gyre |
---|
| 440 | |
---|
| 441 | |
---|
[3] | 442 | SUBROUTINE istate_uvg |
---|
| 443 | !!---------------------------------------------------------------------- |
---|
| 444 | !! *** ROUTINE istate_uvg *** |
---|
| 445 | !! |
---|
| 446 | !! ** Purpose : Compute the geostrophic velocities from (tn,sn) fields |
---|
| 447 | !! |
---|
| 448 | !! ** Method : Using the hydrostatic hypothesis the now hydrostatic |
---|
| 449 | !! pressure is computed by integrating the in-situ density from the |
---|
| 450 | !! surface to the bottom. |
---|
| 451 | !! p=integral [ rau*g dz ] |
---|
| 452 | !!---------------------------------------------------------------------- |
---|
[359] | 453 | USE dynspg ! surface pressure gradient (dyn_spg routine) |
---|
[3] | 454 | USE divcur ! hor. divergence & rel. vorticity (div_cur routine) |
---|
| 455 | USE lbclnk ! ocean lateral boundary condition (or mpp link) |
---|
| 456 | |
---|
| 457 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 458 | INTEGER :: indic ! ??? |
---|
[508] | 459 | REAL(wp) :: zmsv, zphv, zmsu, zphu, zalfg ! temporary scalars |
---|
| 460 | REAL(wp), DIMENSION (jpi,jpj,jpk) :: zprn ! workspace |
---|
[3] | 461 | !!---------------------------------------------------------------------- |
---|
| 462 | |
---|
| 463 | IF(lwp) WRITE(numout,*) |
---|
| 464 | IF(lwp) WRITE(numout,*) 'istate_uvg : Start from Geostrophy' |
---|
| 465 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~' |
---|
| 466 | |
---|
| 467 | ! Compute the now hydrostatic pressure |
---|
| 468 | ! ------------------------------------ |
---|
| 469 | |
---|
[15] | 470 | zalfg = 0.5 * grav * rau0 |
---|
[508] | 471 | |
---|
| 472 | zprn(:,:,1) = zalfg * fse3w(:,:,1) * ( 1 + rhd(:,:,1) ) ! Surface value |
---|
[3] | 473 | |
---|
[508] | 474 | DO jk = 2, jpkm1 ! Vertical integration from the surface |
---|
[3] | 475 | zprn(:,:,jk) = zprn(:,:,jk-1) & |
---|
[359] | 476 | & + zalfg * fse3w(:,:,jk) * ( 2. + rhd(:,:,jk) + rhd(:,:,jk-1) ) |
---|
[3] | 477 | END DO |
---|
| 478 | |
---|
| 479 | ! Compute geostrophic balance |
---|
| 480 | ! --------------------------- |
---|
| 481 | DO jk = 1, jpkm1 |
---|
| 482 | DO jj = 2, jpjm1 |
---|
| 483 | DO ji = fs_2, fs_jpim1 ! vertor opt. |
---|
| 484 | zmsv = 1. / MAX( umask(ji-1,jj+1,jk) + umask(ji ,jj+1,jk) & |
---|
| 485 | + umask(ji-1,jj ,jk) + umask(ji ,jj ,jk) , 1. ) |
---|
| 486 | zphv = ( zprn(ji ,jj+1,jk) - zprn(ji-1,jj+1,jk) ) * umask(ji-1,jj+1,jk) / e1u(ji-1,jj+1) & |
---|
| 487 | + ( zprn(ji+1,jj+1,jk) - zprn(ji ,jj+1,jk) ) * umask(ji ,jj+1,jk) / e1u(ji ,jj+1) & |
---|
| 488 | + ( zprn(ji ,jj ,jk) - zprn(ji-1,jj ,jk) ) * umask(ji-1,jj ,jk) / e1u(ji-1,jj ) & |
---|
| 489 | + ( zprn(ji+1,jj ,jk) - zprn(ji ,jj ,jk) ) * umask(ji ,jj ,jk) / e1u(ji ,jj ) |
---|
| 490 | zphv = 1. / rau0 * zphv * zmsv * vmask(ji,jj,jk) |
---|
| 491 | |
---|
| 492 | zmsu = 1. / MAX( vmask(ji+1,jj ,jk) + vmask(ji ,jj ,jk) & |
---|
| 493 | + vmask(ji+1,jj-1,jk) + vmask(ji ,jj-1,jk) , 1. ) |
---|
| 494 | zphu = ( zprn(ji+1,jj+1,jk) - zprn(ji+1,jj ,jk) ) * vmask(ji+1,jj ,jk) / e2v(ji+1,jj ) & |
---|
| 495 | + ( zprn(ji ,jj+1,jk) - zprn(ji ,jj ,jk) ) * vmask(ji ,jj ,jk) / e2v(ji ,jj ) & |
---|
| 496 | + ( zprn(ji+1,jj ,jk) - zprn(ji+1,jj-1,jk) ) * vmask(ji+1,jj-1,jk) / e2v(ji+1,jj-1) & |
---|
| 497 | + ( zprn(ji ,jj ,jk) - zprn(ji ,jj-1,jk) ) * vmask(ji ,jj-1,jk) / e2v(ji ,jj-1) |
---|
| 498 | zphu = 1. / rau0 * zphu * zmsu * umask(ji,jj,jk) |
---|
| 499 | |
---|
| 500 | ! Compute the geostrophic velocities |
---|
| 501 | un(ji,jj,jk) = -2. * zphu / ( ff(ji,jj) + ff(ji ,jj-1) ) |
---|
| 502 | vn(ji,jj,jk) = 2. * zphv / ( ff(ji,jj) + ff(ji-1,jj ) ) |
---|
| 503 | END DO |
---|
| 504 | END DO |
---|
| 505 | END DO |
---|
| 506 | |
---|
| 507 | IF(lwp) WRITE(numout,*) ' we force to zero bottom velocity' |
---|
| 508 | |
---|
| 509 | ! Susbtract the bottom velocity (level jpk-1 for flat bottom case) |
---|
| 510 | ! to have a zero bottom velocity |
---|
| 511 | |
---|
| 512 | DO jk = 1, jpkm1 |
---|
| 513 | un(:,:,jk) = ( un(:,:,jk) - un(:,:,jpkm1) ) * umask(:,:,jk) |
---|
| 514 | vn(:,:,jk) = ( vn(:,:,jk) - vn(:,:,jpkm1) ) * vmask(:,:,jk) |
---|
| 515 | END DO |
---|
| 516 | |
---|
| 517 | CALL lbc_lnk( un, 'U', -1. ) |
---|
| 518 | CALL lbc_lnk( vn, 'V', -1. ) |
---|
| 519 | |
---|
| 520 | ub(:,:,:) = un(:,:,:) |
---|
| 521 | vb(:,:,:) = vn(:,:,:) |
---|
| 522 | |
---|
| 523 | ! WARNING !!!!! |
---|
| 524 | ! after initializing u and v, we need to calculate the initial streamfunction bsf. |
---|
| 525 | ! Otherwise, only the trend will be computed and the model will blow up (inconsistency). |
---|
| 526 | ! to do that, we call dyn_spg with a special trick: |
---|
[508] | 527 | ! we fill ua and va with the velocities divided by dt, and the streamfunction will be brought to the |
---|
| 528 | ! right value assuming the velocities have been set up in one time step. |
---|
| 529 | ! we then set bsfd to zero (first guess for next step is d(psi)/dt = 0.) |
---|
| 530 | ! sets up s false trend to calculate the barotropic streamfunction. |
---|
[3] | 531 | |
---|
| 532 | ua(:,:,:) = ub(:,:,:) / rdt |
---|
| 533 | va(:,:,:) = vb(:,:,:) / rdt |
---|
| 534 | |
---|
[359] | 535 | ! calls dyn_spg. we assume euler time step, starting from rest. |
---|
[3] | 536 | indic = 0 |
---|
[359] | 537 | CALL dyn_spg( nit000, indic ) ! surface pressure gradient |
---|
[3] | 538 | |
---|
| 539 | ! the new velocity is ua*rdt |
---|
| 540 | |
---|
| 541 | CALL lbc_lnk( ua, 'U', -1. ) |
---|
| 542 | CALL lbc_lnk( va, 'V', -1. ) |
---|
| 543 | |
---|
| 544 | ub(:,:,:) = ua(:,:,:) * rdt |
---|
| 545 | vb(:,:,:) = va(:,:,:) * rdt |
---|
| 546 | ua(:,:,:) = 0.e0 |
---|
| 547 | va(:,:,:) = 0.e0 |
---|
| 548 | un(:,:,:) = ub(:,:,:) |
---|
| 549 | vn(:,:,:) = vb(:,:,:) |
---|
| 550 | |
---|
| 551 | ! Compute the divergence and curl |
---|
| 552 | |
---|
| 553 | CALL div_cur( nit000 ) ! now horizontal divergence and curl |
---|
| 554 | |
---|
| 555 | hdivb(:,:,:) = hdivn(:,:,:) ! set the before to the now value |
---|
| 556 | rotb (:,:,:) = rotn (:,:,:) ! set the before to the now value |
---|
[508] | 557 | ! |
---|
[3] | 558 | END SUBROUTINE istate_uvg |
---|
| 559 | |
---|
| 560 | !!===================================================================== |
---|
| 561 | END MODULE istate |
---|