[325] | 1 | MODULE dtadyn |
---|
| 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE dtadyn *** |
---|
| 4 | !! OFFLINE : interpolation of the physical fields |
---|
| 5 | !!===================================================================== |
---|
| 6 | |
---|
| 7 | !!---------------------------------------------------------------------- |
---|
| 8 | !! dta_dyn_init : initialization, namelist read, and parameters control |
---|
| 9 | !! dta_dyn : Interpolation of the fields |
---|
| 10 | !!---------------------------------------------------------------------- |
---|
| 11 | !! * Modules used |
---|
| 12 | USE oce ! ocean dynamics and tracers variables |
---|
| 13 | USE dom_oce ! ocean space and time domain variables |
---|
| 14 | USE zdf_oce ! ocean vertical physics |
---|
| 15 | USE in_out_manager ! I/O manager |
---|
| 16 | USE phycst ! physical constants |
---|
| 17 | USE ocesbc |
---|
| 18 | USE ldfslp |
---|
| 19 | USE blk_oce |
---|
[446] | 20 | USE ldfeiv ! eddy induced velocity coef. (ldf_eiv routine) |
---|
| 21 | USE ldftra_oce ! ocean tracer lateral physics |
---|
[325] | 22 | USE zdfmxl |
---|
| 23 | USE trabbl ! tracers: bottom boundary layer |
---|
| 24 | USE ocfzpt |
---|
| 25 | USE zdfddm ! vertical physics: double diffusion |
---|
| 26 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
| 27 | USE lib_mpp ! distributed memory computing library |
---|
| 28 | |
---|
| 29 | IMPLICIT NONE |
---|
| 30 | PRIVATE |
---|
| 31 | |
---|
| 32 | !! * Routine accessibility |
---|
| 33 | PUBLIC dta_dyn_init ! called by opa.F90 |
---|
| 34 | PUBLIC dta_dyn ! called by step.F90 |
---|
| 35 | |
---|
| 36 | !! * Module variables |
---|
| 37 | INTEGER , PUBLIC, PARAMETER :: jpflx = 13 |
---|
| 38 | INTEGER , PUBLIC, PARAMETER :: & |
---|
| 39 | jptaux = 1 , & ! indice in flux for taux |
---|
| 40 | jptauy = 2 , & ! indice in flux for tauy |
---|
| 41 | jpwind = 3 , & ! indice in flux for wspd |
---|
| 42 | jpemp = 4 , & ! indice in flux for E-P |
---|
| 43 | jpice = 5 , & ! indice in flux for ice concentration |
---|
| 44 | jpqsr = 6 ! indice in flux for shortwave heat flux |
---|
| 45 | |
---|
| 46 | LOGICAL , PUBLIC :: & |
---|
| 47 | lperdyn = .TRUE. , & ! boolean for periodic fields or not |
---|
| 48 | lfirdyn = .TRUE. ! boolean for the first call or not |
---|
| 49 | |
---|
| 50 | INTEGER , PUBLIC :: & |
---|
| 51 | ndtadyn = 12 , & ! Number of dat in one year |
---|
| 52 | ndtatot = 12 , & ! Number of data in the input field |
---|
| 53 | nsptint = 1 , & ! type of spatial interpolation |
---|
| 54 | nficdyn = 2 ! number of dynamical fields |
---|
| 55 | |
---|
[495] | 56 | INTEGER :: & |
---|
| 57 | ndyn1, ndyn2 , & |
---|
[325] | 58 | nlecoff = 0 , & ! switch for the first read |
---|
| 59 | numfl_t, numfl_u, & |
---|
[495] | 60 | numfl_v, numfl_w |
---|
[325] | 61 | |
---|
| 62 | |
---|
| 63 | REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: & |
---|
| 64 | tdta , & ! temperature at two consecutive times |
---|
| 65 | sdta , & ! salinity at two consecutive times |
---|
| 66 | udta , & ! zonal velocity at two consecutive times |
---|
| 67 | vdta , & ! meridional velocity at two consecutive times |
---|
| 68 | wdta , & ! vertical velocity at two consecutive times |
---|
| 69 | avtdta ! vertical diffusivity coefficient |
---|
| 70 | |
---|
| 71 | #if defined key_ldfslp |
---|
| 72 | REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: & |
---|
| 73 | uslpdta , & ! zonal isopycnal slopes |
---|
| 74 | vslpdta , & ! meridional isopycnal slopes |
---|
| 75 | wslpidta , & ! zonal diapycnal slopes |
---|
| 76 | wslpjdta ! meridional diapycnal slopes |
---|
| 77 | #endif |
---|
| 78 | |
---|
[495] | 79 | #if ! defined key_off_degrad |
---|
| 80 | |
---|
| 81 | # if defined key_traldf_c2d |
---|
[446] | 82 | REAL(wp), DIMENSION(jpi,jpj,2) :: & |
---|
[495] | 83 | ahtwdta ! Lateral diffusivity |
---|
| 84 | # if defined key_trcldf_eiv |
---|
| 85 | REAL(wp), DIMENSION(jpi,jpj,2) :: & |
---|
| 86 | aeiwdta ! G&M coefficient |
---|
| 87 | # endif |
---|
| 88 | # endif |
---|
| 89 | |
---|
| 90 | #else |
---|
| 91 | |
---|
| 92 | REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: & |
---|
| 93 | ahtudta, ahtvdta, ahtwdta ! Lateral diffusivity |
---|
| 94 | # if defined key_trcldf_eiv |
---|
| 95 | REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: & |
---|
| 96 | aeiudta, aeivdta, aeiwdta ! G&M coefficient |
---|
| 97 | # endif |
---|
| 98 | |
---|
[446] | 99 | #endif |
---|
[495] | 100 | # if defined key_diaeiv |
---|
| 101 | !! GM Velocity : to be used latter |
---|
| 102 | REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: & |
---|
| 103 | eivudta, eivvdta, eivwdta |
---|
| 104 | # endif |
---|
[446] | 105 | |
---|
[325] | 106 | REAL(wp), DIMENSION(jpi,jpj,jpflx,2) :: & |
---|
| 107 | flxdta ! auxiliary 2-D forcing fields at two consecutive times |
---|
| 108 | REAL(wp), DIMENSION(jpi,jpj,2) :: & |
---|
| 109 | zmxldta ! mixed layer depth at two consecutive times |
---|
| 110 | |
---|
| 111 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
| 112 | REAL(wp), DIMENSION(jpi,jpj,2) :: & |
---|
| 113 | bblxdta , & ! frequency of bbl in the x direction at 2 consecutive times |
---|
| 114 | bblydta ! frequency of bbl in the y direction at 2 consecutive times |
---|
| 115 | #endif |
---|
| 116 | |
---|
| 117 | !! * Substitutions |
---|
| 118 | # include "domzgr_substitute.h90" |
---|
| 119 | # include "vectopt_loop_substitute.h90" |
---|
[343] | 120 | !!---------------------------------------------------------------------- |
---|
| 121 | !! OPA 9.0 , LOCEAN-IPSL (2005) |
---|
[699] | 122 | !! $Id$ |
---|
[343] | 123 | !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt |
---|
| 124 | !!---------------------------------------------------------------------- |
---|
[325] | 125 | |
---|
| 126 | CONTAINS |
---|
| 127 | |
---|
| 128 | SUBROUTINE dta_dyn_init |
---|
| 129 | !!---------------------------------------------------------------------- |
---|
| 130 | !! *** ROUTINE dta_dyn_init *** |
---|
| 131 | !! |
---|
| 132 | !! ** Purpose : initializations of parameters for the interpolation |
---|
| 133 | !! |
---|
| 134 | !! ** Method : |
---|
| 135 | !! |
---|
| 136 | !! History : |
---|
| 137 | !! ! original : 92-01 (M. Imbard: sub domain) |
---|
| 138 | !! ! 98-04 (L.Bopp MA Foujols: slopes for isopyc.) |
---|
| 139 | !! ! 98-05 (L. Bopp read output of coupled run) |
---|
| 140 | !! ! 05-03 (O. Aumont and A. El Moussaoui) F90 |
---|
| 141 | !!---------------------------------------------------------------------- |
---|
| 142 | !! * Modules used |
---|
| 143 | |
---|
| 144 | !! * Local declarations |
---|
| 145 | |
---|
| 146 | |
---|
| 147 | NAMELIST/nam_offdyn/ ndtadyn, ndtatot, nsptint, & |
---|
| 148 | & nficdyn, lperdyn |
---|
| 149 | !!---------------------------------------------------------------------- |
---|
| 150 | |
---|
| 151 | ! Define the dynamical input parameters |
---|
| 152 | ! ====================================== |
---|
| 153 | |
---|
| 154 | ! Read Namelist nam_offdyn : Lateral physics on tracers |
---|
| 155 | REWIND( numnam ) |
---|
| 156 | READ ( numnam, nam_offdyn ) |
---|
| 157 | |
---|
| 158 | IF(lwp) THEN |
---|
| 159 | WRITE(numout,*) |
---|
| 160 | WRITE(numout,*) 'nam_offdyn : offline dynamical selection' |
---|
| 161 | WRITE(numout,*) '~~~~~~~' |
---|
| 162 | WRITE(numout,*) ' Namelist nam_offdyn : set parameters for the lecture of the dynamical fields' |
---|
| 163 | WRITE(numout,*) |
---|
| 164 | WRITE(numout,*) ' number of elements in the FILE for a year ndtadyn = ' , ndtadyn |
---|
| 165 | WRITE(numout,*) ' total number of elements in the FILE ndtatot = ' , ndtatot |
---|
| 166 | WRITE(numout,*) ' type of interpolation nsptint = ' , nsptint |
---|
| 167 | WRITE(numout,*) ' number of dynamics FILE nficdyn = ' , nficdyn |
---|
| 168 | WRITE(numout,*) ' loop on the same FILE lperdyn = ' , lperdyn |
---|
| 169 | WRITE(numout,*) ' ' |
---|
| 170 | ENDIF |
---|
| 171 | |
---|
| 172 | END SUBROUTINE dta_dyn_init |
---|
| 173 | |
---|
| 174 | SUBROUTINE dta_dyn(kt) |
---|
| 175 | !!---------------------------------------------------------------------- |
---|
| 176 | !! *** ROUTINE dta_dyn *** |
---|
| 177 | !! |
---|
| 178 | !! ** Purpose : Prepares dynamics and physics fields from an |
---|
| 179 | !! OPA9 simulation for an off-line simulation |
---|
| 180 | !! for passive tracer |
---|
| 181 | !! |
---|
| 182 | !! ** Method : calculates the position of DATA to read READ DATA |
---|
| 183 | !! (example month changement) computes slopes IF needed |
---|
| 184 | !! interpolates DATA IF needed |
---|
| 185 | !! |
---|
| 186 | !! ** History : |
---|
| 187 | !! ! original : 92-01 (M. Imbard: sub domain) |
---|
| 188 | !! ! addition : 98-04 (L.Bopp MA Foujols: slopes for isopyc.) |
---|
| 189 | !! ! addition : 98-05 (L. Bopp read output of coupled run) |
---|
| 190 | !! ! addition : 05-03 (O. Aumont and A. El Moussaoui) F90 |
---|
[495] | 191 | !! ! addition : 05-12 (C. Ethe) Adapted for DEGINT |
---|
[325] | 192 | !!---------------------------------------------------------------------- |
---|
| 193 | !! * Modules used |
---|
| 194 | USE eosbn2 |
---|
| 195 | |
---|
| 196 | !! * Arguments |
---|
| 197 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
| 198 | |
---|
| 199 | !! * Local declarations |
---|
| 200 | INTEGER :: iper, iperm1, iswap |
---|
| 201 | |
---|
| 202 | REAL(wp) :: zpdtan, zpdtpe, zdemi, zt |
---|
| 203 | REAL(wp) :: zweigh, zweighm1 |
---|
| 204 | |
---|
| 205 | REAL(wp), DIMENSION(jpi,jpj,jpflx) :: & |
---|
| 206 | flx ! auxiliary field for 2-D surface boundary conditions |
---|
| 207 | |
---|
| 208 | |
---|
| 209 | ! 0. Initialization |
---|
| 210 | ! ----------------- |
---|
[446] | 211 | |
---|
[431] | 212 | IF (lfirdyn) THEN |
---|
[446] | 213 | ! |
---|
| 214 | ! time step MUST BE nint000 |
---|
| 215 | ! |
---|
| 216 | IF (kt.ne.nit000) THEN |
---|
| 217 | IF (lwp) THEN |
---|
| 218 | WRITE (numout,*) ' kt MUST BE EQUAL to nit000. kt=',kt & |
---|
| 219 | ,' nit000=',nit000 |
---|
| 220 | END IF |
---|
| 221 | STOP 'dtadyn' |
---|
| 222 | END if |
---|
| 223 | ! Initialize the parameters of the interpolation |
---|
| 224 | CALL dta_dyn_init |
---|
[431] | 225 | ENDIF |
---|
[325] | 226 | |
---|
[446] | 227 | |
---|
[325] | 228 | zpdtan = raass / rdt |
---|
| 229 | zpdtpe = ((zpdtan / FLOAT (ndtadyn))) |
---|
| 230 | zdemi = zpdtpe * 0.5 |
---|
| 231 | zt = (FLOAT (kt) + zdemi ) / zpdtpe |
---|
| 232 | |
---|
| 233 | zweigh = zt - FLOAT(INT(zt)) |
---|
| 234 | zweighm1 = 1. - zweigh |
---|
| 235 | |
---|
| 236 | IF (lperdyn) THEN |
---|
| 237 | iperm1 = MOD(INT(zt),ndtadyn) |
---|
| 238 | ELSE |
---|
[612] | 239 | iperm1 = MOD(INT(zt),(ndtatot-1)) + 1 |
---|
[325] | 240 | ENDIF |
---|
| 241 | iper = iperm1 + 1 |
---|
| 242 | IF (iperm1 == 0) THEN |
---|
| 243 | IF (lperdyn) THEN |
---|
| 244 | iperm1 = ndtadyn |
---|
| 245 | ELSE |
---|
| 246 | IF (lfirdyn) THEN |
---|
| 247 | IF (lwp) THEN |
---|
| 248 | WRITE (numout,*) ' dynamic file is not periodic ' |
---|
| 249 | WRITE (numout,*) ' with or without interpolation, ' |
---|
| 250 | WRITE (numout,*) ' we take the first value' |
---|
| 251 | WRITE (numout,*) ' for the previous period ' |
---|
| 252 | WRITE (numout,*) ' iperm1 = 0 ' |
---|
| 253 | END IF |
---|
| 254 | END IF |
---|
| 255 | END IF |
---|
| 256 | END IF |
---|
| 257 | |
---|
| 258 | iswap = 0 |
---|
| 259 | |
---|
| 260 | ! 1. First call lfirdyn = true |
---|
| 261 | ! ---------------------------- |
---|
| 262 | |
---|
| 263 | IF (lfirdyn) THEN |
---|
| 264 | ! |
---|
| 265 | ! store the information of the period read |
---|
| 266 | ! |
---|
| 267 | ndyn1 = iperm1 |
---|
| 268 | ndyn2 = iper |
---|
| 269 | |
---|
| 270 | IF (lwp) THEN |
---|
| 271 | WRITE (numout,*) & |
---|
| 272 | ' dynamics DATA READ for the period ndyn1 =',ndyn1, & |
---|
| 273 | & ' and for the period ndyn2 = ',ndyn2 |
---|
| 274 | WRITE (numout,*) ' time step is :',kt |
---|
| 275 | WRITE (numout,*) ' we have ndtadyn = ',ndtadyn,& |
---|
| 276 | & ' records in the dynamic FILE for one year' |
---|
| 277 | END IF |
---|
| 278 | ! |
---|
| 279 | ! DATA READ for the iperm1 period |
---|
| 280 | ! |
---|
[495] | 281 | IF( iperm1 /= 0 ) THEN |
---|
[325] | 282 | CALL dynrea( kt, iperm1 ) |
---|
| 283 | ELSE |
---|
| 284 | CALL dynrea( kt, 1 ) |
---|
| 285 | ENDIF |
---|
| 286 | ! |
---|
| 287 | ! Computes dynamical fields |
---|
| 288 | ! |
---|
| 289 | tn(:,:,:)=tdta(:,:,:,2) |
---|
| 290 | sn(:,:,:)=sdta(:,:,:,2) |
---|
| 291 | avt(:,:,:)=avtdta(:,:,:,2) |
---|
| 292 | |
---|
[495] | 293 | |
---|
[325] | 294 | IF(lwp) THEN |
---|
| 295 | WRITE(numout,*)' temperature ' |
---|
| 296 | WRITE(numout,*) |
---|
[446] | 297 | CALL prihre(tn(1,1,1),jpi,jpj,1,jpi,20,1,jpj,20,1.,numout) |
---|
[325] | 298 | WRITE(numout,*) ' level = ',jpk/2 |
---|
[446] | 299 | CALL prihre(tn(1,1,jpk/2),jpi,jpj,1,jpi,20,1,jpj,20,1.,numout) |
---|
[325] | 300 | WRITE(numout,*) ' level = ',jpkm1 |
---|
[446] | 301 | CALL prihre(tn(1,1,jpkm1),jpi,jpj,1,jpi,20,1,jpj,20,1.,numout) |
---|
[325] | 302 | ENDIF |
---|
| 303 | |
---|
| 304 | #if defined key_ldfslp |
---|
| 305 | CALL eos( tn, sn, rhd, rhop ) ! Time-filtered in situ density |
---|
| 306 | CALL bn2( tn, sn, rn2 ) ! before Brunt-Vaisala frequency |
---|
| 307 | CALL zdf_mxl( kt ) ! mixed layer depth |
---|
| 308 | CALL ldf_slp( kt, rhd, rn2 ) |
---|
| 309 | |
---|
| 310 | uslpdta(:,:,:,2)=uslp(:,:,:) |
---|
| 311 | vslpdta(:,:,:,2)=vslp(:,:,:) |
---|
| 312 | wslpidta(:,:,:,2)=wslpi(:,:,:) |
---|
| 313 | wslpjdta(:,:,:,2)=wslpj(:,:,:) |
---|
| 314 | #endif |
---|
| 315 | ! |
---|
| 316 | ! swap from record 2 to 1 |
---|
| 317 | ! |
---|
| 318 | udta(:,:,:,1)=udta(:,:,:,2) |
---|
| 319 | vdta(:,:,:,1)=vdta(:,:,:,2) |
---|
| 320 | wdta(:,:,:,1)=wdta(:,:,:,2) |
---|
| 321 | avtdta(:,:,:,1)=avtdta(:,:,:,2) |
---|
| 322 | tdta(:,:,:,1)=tdta(:,:,:,2) |
---|
| 323 | sdta(:,:,:,1)=sdta(:,:,:,2) |
---|
| 324 | #if defined key_ldfslp |
---|
| 325 | uslpdta(:,:,:,1)=uslpdta(:,:,:,2) |
---|
| 326 | vslpdta(:,:,:,1)=vslpdta(:,:,:,2) |
---|
| 327 | wslpidta(:,:,:,1)=wslpidta(:,:,:,2) |
---|
| 328 | wslpjdta(:,:,:,1)=wslpjdta(:,:,:,2) |
---|
| 329 | #endif |
---|
| 330 | flxdta(:,:,:,1) = flxdta(:,:,:,2) |
---|
| 331 | zmxldta(:,:,1)=zmxldta(:,:,2) |
---|
[495] | 332 | #if ! defined key_off_degrad |
---|
| 333 | |
---|
| 334 | # if defined key_traldf_c2d |
---|
| 335 | ahtwdta(:,:,1)= ahtwdta(:,:,2) |
---|
| 336 | # if defined key_trcldf_eiv |
---|
| 337 | aeiwdta(:,:,1)= aeiwdta(:,:,2) |
---|
| 338 | # endif |
---|
| 339 | # endif |
---|
| 340 | |
---|
| 341 | #else |
---|
| 342 | ahtudta(:,:,:,1) = ahtudta(:,:,:,2) |
---|
| 343 | ahtvdta(:,:,:,1) = ahtvdta(:,:,:,2) |
---|
| 344 | ahtwdta(:,:,:,1) = ahtwdta(:,:,:,2) |
---|
| 345 | # if defined key_trcldf_eiv |
---|
| 346 | aeiudta(:,:,:,1) = aeiudta(:,:,:,2) |
---|
| 347 | aeivdta(:,:,:,1) = aeivdta(:,:,:,2) |
---|
| 348 | aeiwdta(:,:,:,1) = aeiwdta(:,:,:,2) |
---|
| 349 | # endif |
---|
| 350 | |
---|
[446] | 351 | #endif |
---|
[495] | 352 | |
---|
[325] | 353 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
| 354 | bblxdta(:,:,1)=bblxdta(:,:,2) |
---|
| 355 | bblydta(:,:,1)=bblydta(:,:,2) |
---|
| 356 | #endif |
---|
| 357 | ! |
---|
| 358 | ! indicates a swap |
---|
| 359 | ! |
---|
| 360 | iswap = 1 |
---|
| 361 | ! |
---|
| 362 | ! DATA READ for the iper period |
---|
| 363 | ! |
---|
[495] | 364 | CALL dynrea( kt, iper ) |
---|
[325] | 365 | ! |
---|
| 366 | ! Computes wdta (and slopes if key_trahdfiso) |
---|
| 367 | ! |
---|
| 368 | tn(:,:,:)=tdta(:,:,:,2) |
---|
| 369 | sn(:,:,:)=sdta(:,:,:,2) |
---|
| 370 | avt(:,:,:)=avtdta(:,:,:,2) |
---|
| 371 | |
---|
| 372 | |
---|
| 373 | #if defined key_ldfslp |
---|
| 374 | CALL eos( tn, sn, rhd, rhop ) ! Time-filtered in situ density |
---|
| 375 | CALL bn2( tn, sn, rn2 ) ! before Brunt-Vaisala frequency |
---|
| 376 | CALL zdf_mxl( kt ) ! mixed layer depth |
---|
| 377 | CALL ldf_slp( kt, rhd, rn2 ) |
---|
| 378 | |
---|
| 379 | uslpdta(:,:,:,2)=uslp(:,:,:) |
---|
| 380 | vslpdta(:,:,:,2)=vslp(:,:,:) |
---|
| 381 | wslpidta(:,:,:,2)=wslpi(:,:,:) |
---|
| 382 | wslpjdta(:,:,:,2)=wslpj(:,:,:) |
---|
| 383 | #endif |
---|
| 384 | ! |
---|
| 385 | ! trace the first CALL |
---|
| 386 | ! |
---|
| 387 | lfirdyn=.FALSE. |
---|
| 388 | ENDIF |
---|
| 389 | ! |
---|
| 390 | ! and now what we have to DO at every time step |
---|
| 391 | ! |
---|
| 392 | ! check the validity of the period in memory |
---|
| 393 | ! |
---|
| 394 | IF (iperm1.NE.ndyn1) THEN |
---|
| 395 | IF (iperm1.EQ.0.) THEN |
---|
| 396 | IF (lwp) THEN |
---|
| 397 | WRITE (numout,*) ' dynamic file is not periodic ' |
---|
| 398 | WRITE (numout,*) ' with or without interpolation, ' |
---|
| 399 | WRITE (numout,*) ' we take the last value' |
---|
| 400 | WRITE (numout,*) ' for the last period ' |
---|
| 401 | WRITE (numout,*) ' iperm1 = 12 ' |
---|
| 402 | WRITE (numout,*) ' iper = 13' |
---|
| 403 | ENDIF |
---|
| 404 | iperm1 = 12 |
---|
| 405 | iper =13 |
---|
| 406 | ENDIF |
---|
| 407 | ! |
---|
| 408 | ! we have to prepare a NEW READ of DATA |
---|
| 409 | ! |
---|
| 410 | ! swap from record 2 to 1 |
---|
| 411 | ! |
---|
[495] | 412 | udta(:,:,:,1) = udta(:,:,:,2) |
---|
| 413 | vdta(:,:,:,1) = vdta(:,:,:,2) |
---|
| 414 | wdta(:,:,:,1)= wdta(:,:,:,2) |
---|
| 415 | avtdta(:,:,:,1) = avtdta(:,:,:,2) |
---|
| 416 | tdta(:,:,:,1) = tdta(:,:,:,2) |
---|
| 417 | sdta(:,:,:,1) = sdta(:,:,:,2) |
---|
[325] | 418 | #if defined key_ldfslp |
---|
[495] | 419 | uslpdta(:,:,:,1) = uslpdta(:,:,:,2) |
---|
| 420 | vslpdta(:,:,:,1) = vslpdta(:,:,:,2) |
---|
| 421 | wslpidta(:,:,:,1) = wslpidta(:,:,:,2) |
---|
| 422 | wslpjdta(:,:,:,1) = wslpjdta(:,:,:,2) |
---|
[325] | 423 | #endif |
---|
| 424 | flxdta(:,:,:,1) = flxdta(:,:,:,2) |
---|
[495] | 425 | zmxldta(:,:,1) = zmxldta(:,:,2) |
---|
| 426 | |
---|
| 427 | #if ! defined key_off_degrad |
---|
| 428 | |
---|
| 429 | # if defined key_traldf_c2d |
---|
| 430 | ahtwdta(:,:,1)= ahtwdta(:,:,2) |
---|
| 431 | # if defined key_trcldf_eiv |
---|
| 432 | aeiwdta(:,:,1)= aeiwdta(:,:,2) |
---|
| 433 | # endif |
---|
| 434 | # endif |
---|
| 435 | |
---|
| 436 | #else |
---|
| 437 | ahtudta(:,:,:,1) = ahtudta(:,:,:,2) |
---|
| 438 | ahtvdta(:,:,:,1) = ahtvdta(:,:,:,2) |
---|
| 439 | ahtwdta(:,:,:,1) = ahtwdta(:,:,:,2) |
---|
| 440 | # if defined key_trcldf_eiv |
---|
| 441 | aeiudta(:,:,:,1) = aeiudta(:,:,:,2) |
---|
| 442 | aeivdta(:,:,:,1) = aeivdta(:,:,:,2) |
---|
| 443 | aeiwdta(:,:,:,1) = aeiwdta(:,:,:,2) |
---|
| 444 | # endif |
---|
| 445 | |
---|
[446] | 446 | #endif |
---|
[495] | 447 | |
---|
[325] | 448 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
[495] | 449 | bblxdta(:,:,1) = bblxdta(:,:,2) |
---|
| 450 | bblydta(:,:,1) = bblydta(:,:,2) |
---|
[325] | 451 | #endif |
---|
| 452 | ! |
---|
| 453 | ! indicates a swap |
---|
| 454 | ! |
---|
| 455 | iswap = 1 |
---|
| 456 | ! |
---|
| 457 | ! READ DATA for the iper period |
---|
| 458 | ! |
---|
[495] | 459 | CALL dynrea( kt, iper ) |
---|
[325] | 460 | ! |
---|
| 461 | ! Computes wdta (and slopes if key_trahdfiso) |
---|
| 462 | ! |
---|
| 463 | tn(:,:,:)=tdta(:,:,:,2) |
---|
| 464 | sn(:,:,:)=sdta(:,:,:,2) |
---|
| 465 | avt(:,:,:)=avtdta(:,:,:,2) |
---|
| 466 | |
---|
| 467 | #if defined key_ldfslp |
---|
| 468 | CALL eos( tn, sn, rhd, rhop ) ! Time-filtered in situ density |
---|
| 469 | CALL bn2( tn, sn, rn2 ) ! before Brunt-Vaisala frequency |
---|
| 470 | CALL zdf_mxl( kt ) ! mixed layer depth |
---|
| 471 | CALL ldf_slp( kt, rhd, rn2 ) |
---|
| 472 | |
---|
| 473 | uslpdta(:,:,:,2)=uslp(:,:,:) |
---|
| 474 | vslpdta(:,:,:,2)=vslp(:,:,:) |
---|
| 475 | wslpidta(:,:,:,2)=wslpi(:,:,:) |
---|
| 476 | wslpjdta(:,:,:,2)=wslpj(:,:,:) |
---|
| 477 | #endif |
---|
| 478 | ! |
---|
| 479 | ! store the information of the period read |
---|
| 480 | ! |
---|
| 481 | ndyn1 = ndyn2 |
---|
| 482 | ndyn2 = iper |
---|
| 483 | ! |
---|
[495] | 484 | ! we have READ another period of DATA ! |
---|
[325] | 485 | IF (lwp) THEN |
---|
| 486 | WRITE (numout,*) ' dynamics DATA READ for the period ndyn1 =',ndyn1 |
---|
| 487 | WRITE (numout,*) ' and the period ndyn2 = ',ndyn2 |
---|
| 488 | WRITE (numout,*) ' time step is :',kt |
---|
| 489 | END IF |
---|
| 490 | |
---|
| 491 | END IF |
---|
| 492 | |
---|
| 493 | ! |
---|
| 494 | ! compute the DATA at the given time step |
---|
| 495 | ! |
---|
[495] | 496 | IF ( nsptint == 0 ) THEN |
---|
[325] | 497 | ! |
---|
| 498 | ! no spatial interpolation |
---|
| 499 | ! |
---|
| 500 | ! DATA are probably correct |
---|
| 501 | ! we have to initialize DATA IF we have changed the period |
---|
| 502 | ! |
---|
| 503 | IF (iswap.eq.1) THEN |
---|
| 504 | ! |
---|
| 505 | ! initialize now fields with the NEW DATA READ |
---|
| 506 | ! |
---|
| 507 | un(:,:,:)=udta(:,:,:,2) |
---|
| 508 | vn(:,:,:)=vdta(:,:,:,2) |
---|
| 509 | wn(:,:,:)=wdta(:,:,:,2) |
---|
| 510 | #if defined key_trc_zdfddm |
---|
| 511 | avs(:,:,:)=avtdta(:,:,:,2) |
---|
| 512 | #endif |
---|
| 513 | avt(:,:,:)=avtdta(:,:,:,2) |
---|
| 514 | tn(:,:,:)=tdta(:,:,:,2) |
---|
| 515 | sn(:,:,:)=sdta(:,:,:,2) |
---|
| 516 | #if defined key_ldfslp |
---|
| 517 | uslp(:,:,:)=uslpdta(:,:,:,2) |
---|
| 518 | vslp(:,:,:)=vslpdta(:,:,:,2) |
---|
| 519 | wslpi(:,:,:)=wslpidta(:,:,:,2) |
---|
| 520 | wslpj(:,:,:)=wslpjdta(:,:,:,2) |
---|
| 521 | #endif |
---|
| 522 | flx(:,:,:) = flxdta(:,:,:,2) |
---|
| 523 | hmld(:,:)=zmxldta(:,:,2) |
---|
[495] | 524 | #if ! defined key_off_degrad |
---|
| 525 | |
---|
| 526 | # if defined key_traldf_c2d |
---|
| 527 | ahtwdta(:,:,1)= ahtwdta(:,:,2) |
---|
| 528 | # if defined key_trcldf_eiv |
---|
| 529 | aeiwdta(:,:,1)= aeiwdta(:,:,2) |
---|
| 530 | # endif |
---|
| 531 | # endif |
---|
| 532 | |
---|
| 533 | #else |
---|
| 534 | ahtudta(:,:,:,1) = ahtudta(:,:,:,2) |
---|
| 535 | ahtvdta(:,:,:,1) = ahtvdta(:,:,:,2) |
---|
| 536 | ahtwdta(:,:,:,1) = ahtwdta(:,:,:,2) |
---|
| 537 | # if defined key_trcldf_eiv |
---|
| 538 | aeiudta(:,:,:,1) = aeiudta(:,:,:,2) |
---|
| 539 | aeivdta(:,:,:,1) = aeivdta(:,:,:,2) |
---|
| 540 | aeiwdta(:,:,:,1) = aeiwdta(:,:,:,2) |
---|
| 541 | # endif |
---|
| 542 | |
---|
[446] | 543 | #endif |
---|
[495] | 544 | |
---|
[325] | 545 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
| 546 | bblx(:,:)=bblxdta(:,:,2) |
---|
| 547 | bbly(:,:)=bblydta(:,:,2) |
---|
| 548 | #endif |
---|
| 549 | ! |
---|
| 550 | ! keep needed fluxes |
---|
| 551 | ! |
---|
| 552 | #if defined key_flx_bulk_monthly || defined key_flx_bulk_daily |
---|
| 553 | vatm(:,:) = flx(:,:,jpwind) |
---|
| 554 | #endif |
---|
| 555 | freeze(:,:) = flx(:,:,jpice) |
---|
| 556 | emp(:,:) = flx(:,:,jpemp) |
---|
[446] | 557 | emps(:,:) = emp(:,:) |
---|
[325] | 558 | qsr(:,:) = flx(:,:,jpqsr) |
---|
| 559 | |
---|
| 560 | END IF |
---|
| 561 | |
---|
| 562 | ELSE |
---|
[495] | 563 | IF ( nsptint == 1 ) THEN |
---|
[325] | 564 | ! |
---|
| 565 | ! linear interpolation |
---|
| 566 | ! |
---|
| 567 | ! initialize now fields with a linear interpolation |
---|
| 568 | ! |
---|
| 569 | un(:,:,:) = zweighm1 * udta(:,:,:,1) + zweigh * udta(:,:,:,2) |
---|
| 570 | vn(:,:,:) = zweighm1 * vdta(:,:,:,1) + zweigh * vdta(:,:,:,2) |
---|
| 571 | wn(:,:,:) = zweighm1 * wdta(:,:,:,1) + zweigh * wdta(:,:,:,2) |
---|
| 572 | #if defined key_zdfddm |
---|
| 573 | avs(:,:,:)= zweighm1 * avtdta(:,:,:,1) + zweigh * avtdta(:,:,:,2) |
---|
| 574 | #endif |
---|
| 575 | avt(:,:,:)= zweighm1 * avtdta(:,:,:,1) + zweigh * avtdta(:,:,:,2) |
---|
| 576 | tn(:,:,:) = zweighm1 * tdta(:,:,:,1) + zweigh * tdta(:,:,:,2) |
---|
| 577 | sn(:,:,:) = zweighm1 * sdta(:,:,:,1) + zweigh * sdta(:,:,:,2) |
---|
| 578 | |
---|
| 579 | |
---|
| 580 | #if defined key_ldfslp |
---|
| 581 | uslp(:,:,:) = zweighm1 * uslpdta(:,:,:,1) + zweigh * uslpdta(:,:,:,2) |
---|
| 582 | vslp(:,:,:) = zweighm1 * vslpdta(:,:,:,1) + zweigh * vslpdta(:,:,:,2) |
---|
| 583 | wslpi(:,:,:) = zweighm1 * wslpidta(:,:,:,1) + zweigh * wslpidta(:,:,:,2) |
---|
| 584 | wslpj(:,:,:) = zweighm1 * wslpjdta(:,:,:,1) + zweigh * wslpjdta(:,:,:,2) |
---|
| 585 | #endif |
---|
| 586 | flx(:,:,:) = zweighm1 * flxdta(:,:,:,1) + zweigh * flxdta(:,:,:,2) |
---|
| 587 | hmld(:,:) = zweighm1 * zmxldta(:,:,1) + zweigh * zmxldta(:,:,2) |
---|
[495] | 588 | #if ! defined key_off_degrad |
---|
| 589 | |
---|
| 590 | # if defined key_traldf_c2d |
---|
| 591 | ahtw(:,:) = zweighm1 * ahtwdta(:,:,1) + zweigh * ahtwdta(:,:,2) |
---|
| 592 | # if defined key_trcldf_eiv |
---|
| 593 | aeiw(:,:) = zweighm1 * aeiwdta(:,:,1) + zweigh * aeiwdta(:,:,2) |
---|
| 594 | # endif |
---|
| 595 | # endif |
---|
| 596 | |
---|
| 597 | #else |
---|
| 598 | ahtu(:,:,:) = zweighm1 * ahtudta(:,:,:,1) + zweigh * ahtudta(:,:,:,2) |
---|
| 599 | ahtv(:,:,:) = zweighm1 * ahtvdta(:,:,:,1) + zweigh * ahtvdta(:,:,:,2) |
---|
| 600 | ahtw(:,:,:) = zweighm1 * ahtwdta(:,:,:,1) + zweigh * ahtwdta(:,:,:,2) |
---|
| 601 | # if defined key_trcldf_eiv |
---|
| 602 | aeiu(:,:,:) = zweighm1 * aeiudta(:,:,:,1) + zweigh * aeiudta(:,:,:,2) |
---|
| 603 | aeiv(:,:,:) = zweighm1 * aeivdta(:,:,:,1) + zweigh * aeivdta(:,:,:,2) |
---|
| 604 | aeiw(:,:,:) = zweighm1 * aeiwdta(:,:,:,1) + zweigh * aeiwdta(:,:,:,2) |
---|
| 605 | # endif |
---|
| 606 | |
---|
[446] | 607 | #endif |
---|
[495] | 608 | |
---|
[325] | 609 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
[495] | 610 | bblx(:,:) = zweighm1 * bblxdta(:,:,1) + zweigh * bblxdta(:,:,2) |
---|
| 611 | bbly(:,:) = zweighm1 * bblydta(:,:,1) + zweigh * bblydta(:,:,2) |
---|
[325] | 612 | #endif |
---|
| 613 | ! |
---|
| 614 | ! keep needed fluxes |
---|
| 615 | ! |
---|
| 616 | #if defined key_flx_bulk_monthly || defined key_flx_bulk_daily |
---|
| 617 | vatm(:,:) = flx(:,:,jpwind) |
---|
| 618 | #endif |
---|
| 619 | freeze(:,:) = flx(:,:,jpice) |
---|
[495] | 620 | emp(:,:) = flx(:,:,jpemp) |
---|
| 621 | emps(:,:) = emp(:,:) |
---|
| 622 | qsr(:,:) = flx(:,:,jpqsr) |
---|
[325] | 623 | ! |
---|
| 624 | ! other interpolation |
---|
| 625 | ! |
---|
| 626 | ELSE |
---|
| 627 | |
---|
| 628 | WRITE (numout,*) ' this kind of interpolation don t EXIST' |
---|
| 629 | WRITE (numout,*) ' at the moment. we STOP ' |
---|
| 630 | STOP 'dtadyn' |
---|
| 631 | |
---|
| 632 | END IF |
---|
| 633 | |
---|
| 634 | END IF |
---|
| 635 | ! |
---|
| 636 | ! lb in any case, we need rhop |
---|
| 637 | ! |
---|
| 638 | CALL eos( tn, sn, rhd, rhop ) |
---|
| 639 | |
---|
[495] | 640 | #if ! defined key_off_degrad && defined key_traldf_c2d |
---|
[446] | 641 | ! In case of 2D varying coefficients, we need aeiv and aeiu |
---|
| 642 | IF( lk_traldf_eiv ) CALL ldf_eiv( kt ) ! eddy induced velocity coefficient |
---|
| 643 | #endif |
---|
| 644 | |
---|
[325] | 645 | END SUBROUTINE dta_dyn |
---|
| 646 | |
---|
| 647 | SUBROUTINE dynrea( kt, kenr ) |
---|
| 648 | !!---------------------------------------------------------------------- |
---|
| 649 | !! *** ROUTINE dynrea *** |
---|
| 650 | !! |
---|
| 651 | !! ** Purpose : READ dynamics fiels from OPA9 netcdf output |
---|
| 652 | !! |
---|
| 653 | !! ** Method : READ the kenr records of DATA and store in |
---|
| 654 | !! in udta(...,2), .... |
---|
| 655 | !! |
---|
| 656 | !! ** History : additions : M. Levy et M. Benjelloul jan 2001 |
---|
| 657 | !! (netcdf FORMAT) |
---|
| 658 | !! 05-03 (O. Aumont and A. El Moussaoui) F90 |
---|
[495] | 659 | !! 06-07 : (C. Ethe) use of iom module |
---|
[325] | 660 | !!---------------------------------------------------------------------- |
---|
| 661 | !! * Modules used |
---|
[495] | 662 | USE iom |
---|
[325] | 663 | |
---|
| 664 | !! * Arguments |
---|
| 665 | INTEGER, INTENT( in ) :: kt, kenr ! time index |
---|
| 666 | !! * Local declarations |
---|
[495] | 667 | INTEGER :: ji, jj, jk, jkenr |
---|
[325] | 668 | |
---|
| 669 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: & |
---|
[495] | 670 | zu, zv, zw, zt, zs, zavt , & ! 3-D dynamical fields |
---|
| 671 | zhdiv ! horizontal divergence |
---|
[325] | 672 | |
---|
| 673 | REAL(wp), DIMENSION(jpi,jpj) :: & |
---|
[495] | 674 | zemp, zqsr, zmld, zice, zwspd |
---|
[325] | 675 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
| 676 | REAL(wp), DIMENSION(jpi,jpj) :: & |
---|
| 677 | zbblx, zbbly |
---|
| 678 | #endif |
---|
| 679 | |
---|
[495] | 680 | #if ! defined key_off_degrad |
---|
| 681 | |
---|
| 682 | # if defined key_traldf_c2d |
---|
| 683 | REAL(wp), DIMENSION(jpi,jpj) :: & |
---|
| 684 | zahtw |
---|
| 685 | # if defined key_trcldf_eiv |
---|
| 686 | REAL(wp), DIMENSION(jpi,jpj) :: & |
---|
| 687 | zaeiw |
---|
| 688 | # endif |
---|
| 689 | # endif |
---|
| 690 | |
---|
| 691 | #else |
---|
| 692 | |
---|
| 693 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: & |
---|
| 694 | zahtu, zahtv, zahtw ! Lateral diffusivity |
---|
| 695 | # if defined key_trcldf_eiv |
---|
| 696 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: & |
---|
| 697 | zaeiu, zaeiv, zaeiw ! G&M coefficient |
---|
| 698 | # endif |
---|
| 699 | |
---|
| 700 | #endif |
---|
| 701 | |
---|
| 702 | # if defined key_diaeiv |
---|
| 703 | !! GM Velocity : to be used latter |
---|
| 704 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: & |
---|
| 705 | zeivu, zeivv, zeivw |
---|
| 706 | # endif |
---|
| 707 | |
---|
[325] | 708 | CHARACTER(len=45) :: & |
---|
| 709 | clname_t = 'dyna_grid_T.nc', & |
---|
| 710 | clname_u = 'dyna_grid_U.nc', & |
---|
| 711 | clname_v = 'dyna_grid_V.nc', & |
---|
[495] | 712 | clname_w = 'dyna_grid_W.nc' |
---|
[325] | 713 | ! |
---|
| 714 | ! 0. Initialization |
---|
| 715 | ! cas d'un fichier non periodique : on utilise deux fois le premier et |
---|
| 716 | ! le dernier champ temporel |
---|
| 717 | |
---|
| 718 | jkenr = kenr |
---|
| 719 | |
---|
| 720 | IF(lwp) THEN |
---|
| 721 | WRITE(numout,*) |
---|
| 722 | WRITE(numout,*) 'Dynrea : reading dynamical fields, kenr = ', jkenr |
---|
| 723 | WRITE(numout,*) ' ~~~~~~~' |
---|
[495] | 724 | #if defined key_off_degrad |
---|
| 725 | WRITE(numout,*) ' Degraded fields' |
---|
| 726 | #endif |
---|
[325] | 727 | WRITE(numout,*) |
---|
| 728 | ENDIF |
---|
| 729 | |
---|
| 730 | |
---|
| 731 | IF( kt == nit000 .AND. nlecoff == 0 ) THEN |
---|
| 732 | |
---|
| 733 | nlecoff = 1 |
---|
| 734 | |
---|
[495] | 735 | CALL iom_open ( clname_t, numfl_t ) |
---|
| 736 | CALL iom_open ( clname_u, numfl_u ) |
---|
| 737 | CALL iom_open ( clname_v, numfl_v ) |
---|
| 738 | CALL iom_open ( clname_w, numfl_w ) |
---|
[325] | 739 | |
---|
| 740 | ENDIF |
---|
| 741 | |
---|
[495] | 742 | ! file grid-T |
---|
| 743 | !--------------- |
---|
| 744 | CALL iom_get ( numfl_t, jpdom_data, 'votemper', zt (:,:,:), jkenr ) |
---|
| 745 | CALL iom_get ( numfl_t, jpdom_data, 'vosaline', zs (:,:,:), jkenr ) |
---|
| 746 | CALL iom_get ( numfl_t, jpdom_data, 'somixhgt', zmld (:,: ), jkenr ) |
---|
| 747 | CALL iom_get ( numfl_t, jpdom_data, 'sowaflup', zemp (:,: ), jkenr ) |
---|
| 748 | CALL iom_get ( numfl_t, jpdom_data, 'soshfldo', zqsr (:,: ), jkenr ) |
---|
| 749 | CALL iom_get ( numfl_t, jpdom_data, 'soicecov', zice (:,: ), jkenr ) |
---|
| 750 | CALL iom_get ( numfl_t, jpdom_data, 'sowindsp', zwspd(:,: ), jkenr ) |
---|
[325] | 751 | |
---|
[495] | 752 | ! file grid-U |
---|
| 753 | !--------------- |
---|
| 754 | CALL iom_get ( numfl_u, jpdom_data, 'vozocrtx', zu (:,:,:), jkenr ) |
---|
[325] | 755 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
[495] | 756 | CALL iom_get ( numfl_u, jpdom_data, 'sobblcox', zbblx(:,: ), jkenr ) |
---|
[325] | 757 | #endif |
---|
| 758 | |
---|
[495] | 759 | #if defined key_diaeiv |
---|
| 760 | !! GM Velocity : to be used latter |
---|
| 761 | CALL iom_get ( numfl_u, jpdom_data, 'vozoeivu', zeivu(:,:,:), jkenr ) |
---|
[446] | 762 | #endif |
---|
[325] | 763 | |
---|
[495] | 764 | # if defined key_off_degrad |
---|
| 765 | CALL iom_get ( numfl_u, jpdom_data, 'vozoahtu', zahtu(:,:,:), jkenr ) |
---|
| 766 | # if defined key_trcldf_eiv |
---|
| 767 | CALL iom_get ( numfl_u, jpdom_data, 'vozoaeiu', zaeiu(:,:,:), jkenr ) |
---|
| 768 | # endif |
---|
| 769 | #endif |
---|
[325] | 770 | |
---|
[495] | 771 | ! file grid-V |
---|
| 772 | !--------------- |
---|
| 773 | CALL iom_get ( numfl_v, jpdom_data, 'vomecrty', zv (:,:,:), jkenr ) |
---|
[325] | 774 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
[495] | 775 | CALL iom_get ( numfl_v, jpdom_data, 'sobblcoy', zbbly(:,: ), jkenr ) |
---|
[325] | 776 | #endif |
---|
| 777 | |
---|
[495] | 778 | #if defined key_diaeiv |
---|
| 779 | !! GM Velocity : to be used latter |
---|
| 780 | CALL iom_get ( numfl_v, jpdom_data, 'vomeeivv', zeivv(:,:,:), jkenr ) |
---|
[446] | 781 | #endif |
---|
| 782 | |
---|
[495] | 783 | #if defined key_off_degrad |
---|
| 784 | CALL iom_get ( numfl_v, jpdom_data, 'vomeahtv', zahtv(:,:,:), jkenr ) |
---|
| 785 | # if defined key_trcldf_eiv |
---|
| 786 | CALL iom_get ( numfl_v, jpdom_data, 'vomeaeiv', zaeiv(:,:,:), jkenr ) |
---|
| 787 | # endif |
---|
[446] | 788 | #endif |
---|
[325] | 789 | |
---|
[495] | 790 | ! file grid-W |
---|
| 791 | !--------------- |
---|
| 792 | !! CALL iom_get ( numfl_w, jpdom_data, 'vovecrtz', zw (:,:,:), jkenr ) |
---|
| 793 | # if defined key_zdfddm |
---|
| 794 | CALL iom_get ( numfl_w, jpdom_data, 'voddmavs', zavt (:,:,:), jkenr ) |
---|
[325] | 795 | #else |
---|
[495] | 796 | CALL iom_get ( numfl_w, jpdom_data, 'votkeavt', zavt (:,:,:), jkenr ) |
---|
| 797 | #endif |
---|
[325] | 798 | |
---|
[495] | 799 | # if defined key_diaeiv |
---|
| 800 | !! GM Velocity : to be used latter |
---|
| 801 | CALL iom_get ( numfl_w, jpdom_data, 'voveeivw', zeivw(:,:,:), jkenr ) |
---|
| 802 | #endif |
---|
[446] | 803 | |
---|
[495] | 804 | #if ! defined key_off_degrad |
---|
| 805 | # if defined key_traldf_c2d |
---|
| 806 | CALL iom_get ( numfl_w, jpdom_data, 'soleahtw', zahtw (:,: ), jkenr ) |
---|
[612] | 807 | # if defined key_trcldf_eiv |
---|
[495] | 808 | CALL iom_get ( numfl_w, jpdom_data, 'soleaeiw', zaeiw (:,: ), jkenr ) |
---|
| 809 | # endif |
---|
| 810 | # endif |
---|
| 811 | #else |
---|
| 812 | !! degradation-integration |
---|
| 813 | CALL iom_get ( numfl_w, jpdom_data, 'voveahtw', zahtw(:,:,:), jkenr ) |
---|
| 814 | # if defined key_trcldf_eiv |
---|
| 815 | CALL iom_get ( numfl_w, jpdom_data, 'voveaeiw', zaeiw(:,:,:), jkenr ) |
---|
| 816 | # endif |
---|
[446] | 817 | #endif |
---|
| 818 | |
---|
[495] | 819 | udta(:,:,:,2) = zu(:,:,:) * umask(:,:,:) |
---|
| 820 | vdta(:,:,:,2) = zv(:,:,:) * vmask(:,:,:) |
---|
| 821 | !! wdta(:,:,:,2) = zw(:,:,:) * tmask(:,:,:) |
---|
[325] | 822 | |
---|
| 823 | |
---|
[495] | 824 | ! Computation of vertical velocity using horizontal divergence |
---|
| 825 | zhdiv(:,:,:) = 0. |
---|
| 826 | DO jk = 1, jpkm1 |
---|
| 827 | DO jj = 2, jpjm1 |
---|
| 828 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[612] | 829 | #if defined key_zco |
---|
| 830 | zhdiv(ji,jj,jk) = ( e2u(ji,jj) * udta(ji,jj,jk,2) - e2u(ji-1,jj ) * udta(ji-1,jj ,jk,2) & |
---|
| 831 | & + e1v(ji,jj) * vdta(ji,jj,jk,2) - e1v(ji ,jj-1) * vdta(ji ,jj-1,jk,2) ) & |
---|
[495] | 832 | & / ( e1t(ji,jj) * e2t(ji,jj) ) |
---|
[612] | 833 | #else |
---|
| 834 | zhdiv(ji,jj,jk) = & |
---|
| 835 | ( e2u(ji,jj)*fse3u(ji,jj,jk)*udta(ji,jj,jk,2) - e2u(ji-1,jj)*fse3u(ji-1,jj,jk)*udta(ji-1,jj,jk,2) & |
---|
| 836 | + e1v(ji,jj)*fse3v(ji,jj,jk)*vdta(ji,jj,jk,2) - e1v(ji,jj-1)*fse3v(ji,jj-1,jk)*vdta(ji,jj-1,jk,2) ) & |
---|
| 837 | / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) |
---|
| 838 | #endif |
---|
[495] | 839 | END DO |
---|
[612] | 840 | END DO |
---|
| 841 | ENDDO |
---|
| 842 | |
---|
| 843 | ! Lateral boundary conditions on hdiv |
---|
| 844 | CALL lbc_lnk( zhdiv, 'T', 1. ) |
---|
| 845 | |
---|
| 846 | |
---|
| 847 | ! computation of vertical velocity from the bottom |
---|
[495] | 848 | zw(:,:,jpk) = 0. |
---|
| 849 | DO jk = jpkm1, 1, -1 |
---|
| 850 | zw(:,:,jk) = zw(:,:,jk+1) - fse3t(:,:,jk) * zhdiv(:,:,jk) |
---|
| 851 | END DO |
---|
| 852 | wdta(:,:,:,2) = zw(:,:,:) * tmask(:,:,:) |
---|
[446] | 853 | |
---|
[325] | 854 | |
---|
[495] | 855 | tdta(:,:,:,2) = zt(:,:,:) * tmask(:,:,:) |
---|
| 856 | sdta(:,:,:,2) = zs(:,:,:) * tmask(:,:,:) |
---|
| 857 | avtdta(:,:,:,2) = zavt(:,:,:) * tmask(:,:,:) |
---|
| 858 | #if ! defined key_off_degrad && defined key_traldf_c2d |
---|
| 859 | ahtwdta(:,:,2) = zahtw(:,:) * tmask(:,:,1) |
---|
[612] | 860 | #if defined key_trcldf_eiv |
---|
[495] | 861 | aeiwdta(:,:,2) = zaeiw(:,:) * tmask(:,:,1) |
---|
[446] | 862 | #endif |
---|
| 863 | #endif |
---|
[495] | 864 | |
---|
| 865 | #if defined key_off_degrad |
---|
| 866 | ahtudta(:,:,:,2) = zahtu(:,:,:) * umask(:,:,:) |
---|
| 867 | ahtvdta(:,:,:,2) = zahtv(:,:,:) * vmask(:,:,:) |
---|
| 868 | ahtwdta(:,:,:,2) = zahtw(:,:,:) * tmask(:,:,:) |
---|
| 869 | # if defined key_trcldf_eiv |
---|
| 870 | aeiudta(:,:,:,2) = zaeiu(:,:,:) * umask(:,:,:) |
---|
| 871 | aeivdta(:,:,:,2) = zaeiv(:,:,:) * vmask(:,:,:) |
---|
| 872 | aeiwdta(:,:,:,2) = zaeiw(:,:,:) * tmask(:,:,:) |
---|
| 873 | # endif |
---|
[325] | 874 | #endif |
---|
| 875 | |
---|
| 876 | ! |
---|
| 877 | ! flux : |
---|
| 878 | ! |
---|
[495] | 879 | flxdta(:,:,jpwind,2) = zwspd(:,:) * tmask(:,:,1) |
---|
| 880 | flxdta(:,:,jpice,2) = MIN( 1., zice(:,:) ) * tmask(:,:,1) |
---|
| 881 | flxdta(:,:,jpemp,2) = zemp(:,:) * tmask(:,:,1) |
---|
| 882 | flxdta(:,:,jpqsr,2) = zqsr(:,:) * tmask(:,:,1) |
---|
| 883 | zmxldta(:,:,2) = zmld(:,:) * tmask(:,:,1) |
---|
| 884 | |
---|
[325] | 885 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
[495] | 886 | bblxdta(:,:,2) = MAX( 0., zbblx(:,:) ) |
---|
| 887 | bblydta(:,:,2) = MAX( 0., zbbly(:,:) ) |
---|
[325] | 888 | |
---|
[495] | 889 | WHERE( bblxdta(:,:,2) > 2. ) bblxdta(:,:,2) = 0. |
---|
| 890 | WHERE( bblydta(:,:,2) > 2. ) bblydta(:,:,2) = 0. |
---|
| 891 | |
---|
[325] | 892 | #endif |
---|
| 893 | |
---|
[495] | 894 | IF( kt == nitend ) THEN |
---|
| 895 | CALL iom_close ( numfl_t ) |
---|
| 896 | CALL iom_close ( numfl_u ) |
---|
| 897 | CALL iom_close ( numfl_v ) |
---|
| 898 | CALL iom_close ( numfl_w ) |
---|
| 899 | ENDIF |
---|
| 900 | |
---|
[325] | 901 | END SUBROUTINE dynrea |
---|
| 902 | |
---|
| 903 | |
---|
| 904 | |
---|
| 905 | END MODULE dtadyn |
---|