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