Changeset 1501
- Timestamp:
- 2009-07-20T17:14:38+02:00 (15 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/NEMO/OFF_SRC/dtadyn.F90
r1323 r1501 21 21 USE zdfmxl 22 22 USE trabbl ! tracers: bottom boundary layer 23 USE eosbn2 23 24 USE zdfddm ! vertical physics: double diffusion 24 25 USE lbclnk ! ocean lateral boundary conditions (or mpp link) … … 33 34 PUBLIC dta_dyn ! called by step.F90 34 35 35 !! * Module variables36 INTEGER , PUBLIC, PARAMETER :: jpflx = 737 INTEGER , PUBLIC, PARAMETER :: &38 jptaux = 1 , & ! indice in flux for taux39 jptauy = 2 , & ! indice in flux for tauy40 jpwind = 3 , & ! indice in flux for wind speed41 jpemp = 4 , & ! indice in flux for E-P42 jpice = 5 , & ! indice in flux for ice concentration43 jpqsr = 6 ! indice in flux for shortwave heat flux44 45 36 LOGICAL , PUBLIC :: & 46 37 lperdyn = .TRUE. , & ! boolean for periodic fields or not … … 48 39 49 40 INTEGER , PUBLIC :: & 50 ndtadyn = 12, & ! Number of dat in one year51 ndtatot = 12, & ! Number of data in the input field41 ndtadyn = 73 , & ! Number of dat in one year 42 ndtatot = 73 , & ! Number of data in the input field 52 43 nsptint = 1 , & ! type of spatial interpolation 53 44 nficdyn = 2 ! number of dynamical fields … … 58 49 numfl_t, numfl_u, & 59 50 numfl_v, numfl_w 60 51 52 CHARACTER(len=45) :: & 53 cfile_grid_T = 'dyna_grid_T.nc', & !: name of the grid_T file 54 cfile_grid_U = 'dyna_grid_U.nc', & !: name of the grid_U file 55 cfile_grid_V = 'dyna_grid_V.nc', & !: name of the grid_V file 56 cfile_grid_W = 'dyna_grid_W.nc' !: name of the grid_W file 61 57 62 58 REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: & … … 66 62 vdta , & ! meridional velocity at two consecutive times 67 63 wdta , & ! vertical velocity at two consecutive times 64 #if defined key_trc_diatrd 68 65 hdivdta, & ! horizontal divergence 66 #endif 69 67 avtdta ! vertical diffusivity coefficient 70 68 69 REAL(wp), DIMENSION(jpi,jpj,2) :: & 70 hmlddta, & ! mixed layer depth at two consecutive times 71 wspddta, & ! wind speed at two consecutive times 72 frlddta, & ! sea-ice fraction at two consecutive times 73 empdta , & ! E-P at two consecutive times 74 qsrdta ! short wave heat flux at two consecutive times 71 75 72 76 #if defined key_ldfslp … … 78 82 #endif 79 83 80 #if ! defined key_off_degrad 81 82 # if defined key_traldf_c2d 84 #if ! defined key_off_degrad && defined key_traldf_c2d 83 85 REAL(wp), DIMENSION(jpi,jpj,2) :: & 84 86 ahtwdta ! Lateral diffusivity … … 87 89 aeiwdta ! G&M coefficient 88 90 # endif 89 # endif 90 91 #else 92 91 #endif 92 93 #if defined key_off_degrad 93 94 REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: & 94 95 ahtudta, ahtvdta, ahtwdta ! Lateral diffusivity … … 99 100 100 101 #endif 101 # if defined key_diaeiv102 !! GM Velocity : to be used latter103 REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: &104 eivudta, eivvdta, eivwdta105 # endif106 107 REAL(wp), DIMENSION(jpi,jpj,jpflx,2) :: &108 flxdta ! auxiliary 2-D forcing fields at two consecutive times109 REAL(wp), DIMENSION(jpi,jpj,2) :: &110 zmxldta ! mixed layer depth at two consecutive times111 102 112 103 #if defined key_trcbbl_dif || defined key_trcbbl_adv … … 127 118 CONTAINS 128 119 129 SUBROUTINE dta_dyn_init 130 !!---------------------------------------------------------------------- 131 !! *** ROUTINE dta_dyn_init *** 132 !! 133 !! ** Purpose : initializations of parameters for the interpolation 134 !! 135 !! ** Method : 136 !! 137 !! History : 138 !! ! original : 92-01 (M. Imbard: sub domain) 139 !! ! 98-04 (L.Bopp MA Foujols: slopes for isopyc.) 140 !! ! 98-05 (L. Bopp read output of coupled run) 141 !! ! 05-03 (O. Aumont and A. El Moussaoui) F90 142 !!---------------------------------------------------------------------- 143 !! * Modules used 144 145 !! * Local declarations 146 NAMELIST/namdyn/ ndtadyn, ndtatot, nsptint, & 147 & nficdyn, lperdyn 148 !!---------------------------------------------------------------------- 149 150 ! Define the dynamical input parameters 151 ! ====================================== 152 153 ! Read Namelist namdyn : Lateral physics on tracers 154 REWIND( numnam ) 155 READ ( numnam, namdyn ) 156 157 IF(lwp) THEN 158 WRITE(numout,*) 159 WRITE(numout,*) 'namdyn : offline dynamical selection' 160 WRITE(numout,*) '~~~~~~~' 161 WRITE(numout,*) ' Namelist namdyn : set parameters for the lecture of the dynamical fields' 162 WRITE(numout,*) 163 WRITE(numout,*) ' number of elements in the FILE for a year ndtadyn = ' , ndtadyn 164 WRITE(numout,*) ' total number of elements in the FILE ndtatot = ' , ndtatot 165 WRITE(numout,*) ' type of interpolation nsptint = ' , nsptint 166 WRITE(numout,*) ' number of dynamics FILE nficdyn = ' , nficdyn 167 WRITE(numout,*) ' loop on the same FILE lperdyn = ' , lperdyn 168 WRITE(numout,*) ' ' 169 ENDIF 170 171 END SUBROUTINE dta_dyn_init 172 173 SUBROUTINE dta_dyn(kt) 120 SUBROUTINE dta_dyn( kt ) 174 121 !!---------------------------------------------------------------------- 175 122 !! *** ROUTINE dta_dyn *** … … 190 137 !! ! addition : 05-12 (C. Ethe) Adapted for DEGINT 191 138 !!---------------------------------------------------------------------- 192 !! * Modules used193 USE eosbn2194 195 139 !! * Arguments 196 140 INTEGER, INTENT( in ) :: kt ! ocean time-step index 197 141 198 142 !! * Local declarations 199 INTEGER :: ji,jj,iper, iperm1, iswap143 INTEGER :: iper, iperm1, iswap, izt 200 144 201 145 REAL(wp) :: zpdtan, zpdtpe, zdemi, zt 202 REAL(wp) :: zweigh, zweighm1 203 204 REAL(wp), DIMENSION(jpi,jpj,jpflx) :: & 205 flx ! auxiliary field for 2-D surface boundary conditions 146 REAL(wp) :: zweigh 206 147 207 148 ! 0. Initialization 208 149 ! ----------------- 209 150 210 IF (lfirdyn) THEN151 IF( lfirdyn ) THEN 211 152 ! first time step MUST BE nit000 212 153 IF( kt /= nit000 ) THEN 213 154 IF (lwp) THEN 214 WRITE (numout,*) ' kt MUST BE EQUAL to nit000. kt =',kt ,' nit000=',nit000155 WRITE (numout,*) ' kt MUST BE EQUAL to nit000. kt = ',kt ,' nit000 = ',nit000 215 156 STOP 'dtadyn' 216 157 ENDIF … … 222 163 223 164 zpdtan = raass / rdt 224 zpdtpe = ((zpdtan / FLOAT (ndtadyn)))165 zpdtpe = zpdtan / FLOAT( ndtadyn ) 225 166 zdemi = zpdtpe * 0.5 226 zt = ( FLOAT (kt) + zdemi ) / zpdtpe227 228 zweigh = zt - FLOAT(INT(zt))229 zweigh m1 = 1. - zweigh167 zt = ( FLOAT (kt) + zdemi ) / zpdtpe 168 169 izt = INT( zt ) 170 zweigh = zt - FLOAT( INT(zt) ) 230 171 231 172 IF( lperdyn ) THEN 232 iperm1 = MOD( INT(zt),ndtadyn)173 iperm1 = MOD( izt, ndtadyn ) 233 174 ELSE 234 iperm1 = MOD( INT(zt),(ndtatot-1)) + 1175 iperm1 = MOD( izt, ndtatot - 1 ) + 1 235 176 ENDIF 177 236 178 iper = iperm1 + 1 237 179 IF( iperm1 == 0 ) THEN … … 252 194 ! ---------------------------- 253 195 254 IF (lfirdyn) THEN 255 ! 256 ! store the information of the period read 257 ! 258 ndyn1 = iperm1 259 ndyn2 = iper 260 261 IF (lwp) THEN 262 WRITE (numout,*) & 263 ' dynamics DATA READ for the period ndyn1 =',ndyn1, & 264 & ' and for the period ndyn2 = ',ndyn2 265 WRITE (numout,*) ' time step is :',kt 266 WRITE (numout,*) ' we have ndtadyn = ',ndtadyn,& 267 & ' records in the dynamic FILE for one year' 268 END IF 269 ! 270 ! DATA READ for the iperm1 period 271 ! 272 IF( iperm1 /= 0 ) THEN 273 CALL dynrea( kt, iperm1 ) 274 ELSE 275 CALL dynrea( kt, 1 ) 276 ENDIF 277 ! 278 ! Computes dynamical fields 279 ! 280 tn(:,:,:)=tdta(:,:,:,2) 281 sn(:,:,:)=sdta(:,:,:,2) 282 avt(:,:,:)=avtdta(:,:,:,2) 283 284 285 IF(lwp) THEN 286 WRITE(numout,*)' temperature ' 287 WRITE(numout,*) 288 CALL prihre(tn(1,1,1),jpi,jpj,1,jpi,20,1,jpj,20,1.,numout) 289 WRITE(numout,*) ' level = ',jpk/2 290 CALL prihre(tn(1,1,jpk/2),jpi,jpj,1,jpi,20,1,jpj,20,1.,numout) 291 WRITE(numout,*) ' level = ',jpkm1 292 CALL prihre(tn(1,1,jpkm1),jpi,jpj,1,jpi,20,1,jpj,20,1.,numout) 293 ENDIF 294 196 IF( lfirdyn ) THEN 197 ! store the information of the period read 198 ndyn1 = iperm1 199 ndyn2 = iper 200 201 IF (lwp) THEN 202 WRITE (numout,*) ' dynamics data read for the period ndyn1 =',ndyn1, & 203 & ' and for the period ndyn2 = ',ndyn2 204 WRITE (numout,*) ' time step is : ', kt 205 WRITE (numout,*) ' we have ndtadyn = ',ndtadyn,' records in the dynamic file for one year' 206 END IF 207 ! 208 IF( iperm1 /= 0 ) THEN ! data read for the iperm1 period 209 CALL dynrea( kt, iperm1 ) 210 ELSE 211 CALL dynrea( kt, 1 ) 212 ENDIF 213 295 214 #if defined key_ldfslp 296 CALL eos( tn, sn, rhd, rhop ) ! Time-filtered in situ density 297 CALL bn2( tn, sn, rn2 ) ! before Brunt-Vaisala frequency 298 IF( ln_zps ) & 299 CALL zps_hde( kt, tn , sn , rhd, & ! Partial steps: before Horizontal DErivative 300 gtu, gsu, gru, & ! of t, s, rd at the bottom ocean level 301 gtv, gsv, grv ) 302 CALL zdf_mxl( kt ) ! mixed layer depth 303 CALL ldf_slp( kt, rhd, rn2 ) 304 305 uslpdta (:,:,:,2) = uslp(:,:,:) 306 vslpdta (:,:,:,2) = vslp(:,:,:) 307 wslpidta(:,:,:,2) = wslpi(:,:,:) 308 wslpjdta(:,:,:,2) = wslpj(:,:,:) 309 #endif 310 ! 311 ! swap from record 2 to 1 312 ! 313 udta(:,:,:,1)=udta(:,:,:,2) 314 vdta(:,:,:,1)=vdta(:,:,:,2) 315 wdta(:,:,:,1)=wdta(:,:,:,2) 316 #if defined key_trc_diatrd 317 hdivdta(:,:,:,1)=hdivdta(:,:,:,2) 318 #endif 319 avtdta(:,:,:,1)=avtdta(:,:,:,2) 320 tdta(:,:,:,1)=tdta(:,:,:,2) 321 sdta(:,:,:,1)=sdta(:,:,:,2) 215 ! Computes slopes 216 ! Caution : here tn, sn and avt are used as workspace 217 tn (:,:,:) = tdta (:,:,:,2) 218 sn (:,:,:) = sdta (:,:,:,2) 219 avt(:,:,:) = avtdta(:,:,:,2) 220 221 CALL eos( tn, sn, rhd, rhop ) ! Time-filtered in situ density 222 CALL bn2( tn, sn, rn2 ) ! before Brunt-Vaisala frequency 223 IF( ln_zps ) & 224 & CALL zps_hde( kt, tn , sn , rhd, & ! Partial steps: before Horizontal DErivative 225 & gtu, gsu, gru, & ! of t, s, rd at the bottom ocean level 226 & gtv, gsv, grv ) 227 CALL zdf_mxl( kt ) ! mixed layer depth 228 CALL ldf_slp( kt, rhd, rn2 ) 229 230 uslpdta (:,:,:,2) = uslp (:,:,:) 231 vslpdta (:,:,:,2) = vslp (:,:,:) 232 wslpidta(:,:,:,2) = wslpi(:,:,:) 233 wslpjdta(:,:,:,2) = wslpj(:,:,:) 234 #endif 235 236 ! swap from record 2 to 1 237 CALL swap_dyn_data 238 239 iswap = 1 ! indicates swap 240 241 CALL dynrea( kt, iper ) ! data read for the iper period 242 322 243 #if defined key_ldfslp 323 uslpdta(:,:,:,1)=uslpdta(:,:,:,2) 324 vslpdta(:,:,:,1)=vslpdta(:,:,:,2) 325 wslpidta(:,:,:,1)=wslpidta(:,:,:,2) 326 wslpjdta(:,:,:,1)=wslpjdta(:,:,:,2) 327 #endif 328 flxdta(:,:,:,1) = flxdta(:,:,:,2) 329 zmxldta(:,:,1)=zmxldta(:,:,2) 330 #if ! defined key_off_degrad 331 332 # if defined key_traldf_c2d 333 ahtwdta(:,:,1)= ahtwdta(:,:,2) 334 # if defined key_trcldf_eiv 335 aeiwdta(:,:,1)= aeiwdta(:,:,2) 336 # endif 337 # endif 338 339 #else 340 ahtudta(:,:,:,1) = ahtudta(:,:,:,2) 341 ahtvdta(:,:,:,1) = ahtvdta(:,:,:,2) 342 ahtwdta(:,:,:,1) = ahtwdta(:,:,:,2) 343 # if defined key_trcldf_eiv 344 aeiudta(:,:,:,1) = aeiudta(:,:,:,2) 345 aeivdta(:,:,:,1) = aeivdta(:,:,:,2) 346 aeiwdta(:,:,:,1) = aeiwdta(:,:,:,2) 347 # endif 348 349 #endif 350 351 #if defined key_trcbbl_dif || defined key_trcbbl_adv 352 bblxdta(:,:,1)=bblxdta(:,:,2) 353 bblydta(:,:,1)=bblydta(:,:,2) 354 #endif 355 ! 356 ! indicates a swap 357 ! 358 iswap = 1 359 ! 360 ! DATA READ for the iper period 361 ! 362 CALL dynrea( kt, iper ) 363 ! 364 ! Computes wdta (and slopes if key_trahdfiso) 365 ! 366 tn(:,:,:)=tdta(:,:,:,2) 367 sn(:,:,:)=sdta(:,:,:,2) 368 avt(:,:,:)=avtdta(:,:,:,2) 369 370 371 #if defined key_ldfslp 372 CALL eos( tn, sn, rhd, rhop ) ! Time-filtered in situ density 373 CALL bn2( tn, sn, rn2 ) ! before Brunt-Vaisala frequency 374 IF( ln_zps ) & 375 CALL zps_hde( kt, tn , sn , rhd, & ! Partial steps: before Horizontal DErivative 376 gtu, gsu, gru, & ! of t, s, rd at the bottom ocean level 377 gtv, gsv, grv ) 378 CALL zdf_mxl( kt ) ! mixed layer depth 379 CALL ldf_slp( kt, rhd, rn2 ) 380 381 uslpdta(:,:,:,2)=uslp(:,:,:) 382 vslpdta(:,:,:,2)=vslp(:,:,:) 383 wslpidta(:,:,:,2)=wslpi(:,:,:) 384 wslpjdta(:,:,:,2)=wslpj(:,:,:) 385 #endif 386 ! 387 ! trace the first CALL 388 ! 389 lfirdyn=.FALSE. 244 ! Computes slopes 245 ! Caution : here tn, sn and avt are used as workspace 246 tn (:,:,:) = tdta (:,:,:,2) 247 sn (:,:,:) = sdta (:,:,:,2) 248 avt(:,:,:) = avtdta(:,:,:,2) 249 250 CALL eos( tn, sn, rhd, rhop ) ! Time-filtered in situ density 251 CALL bn2( tn, sn, rn2 ) ! before Brunt-Vaisala frequency 252 IF( ln_zps ) & 253 & CALL zps_hde( kt, tn , sn , rhd, & ! Partial steps: before Horizontal DErivative 254 & gtu, gsu, gru, & ! of t, s, rd at the bottom ocean level 255 & gtv, gsv, grv ) 256 CALL zdf_mxl( kt ) ! mixed layer depth 257 CALL ldf_slp( kt, rhd, rn2 ) 258 259 uslpdta (:,:,:,2) = uslp (:,:,:) 260 vslpdta (:,:,:,2) = vslp (:,:,:) 261 wslpidta(:,:,:,2) = wslpi(:,:,:) 262 wslpjdta(:,:,:,2) = wslpj(:,:,:) 263 #endif 264 ! 265 lfirdyn=.FALSE. ! trace the first call 390 266 ENDIF 391 267 ! 392 ! and now what we have to DO at every time step 393 ! 268 ! And now what we have to do at every time step 394 269 ! check the validity of the period in memory 395 270 ! 396 IF (iperm1.NE.ndyn1) THEN 397 IF (iperm1.EQ.0.) THEN 398 IF (lwp) THEN 399 WRITE (numout,*) ' dynamic file is not periodic ' 400 WRITE (numout,*) ' with or without interpolation, ' 401 WRITE (numout,*) ' we take the last value' 402 WRITE (numout,*) ' for the last period ' 403 WRITE (numout,*) ' iperm1 = 12 ' 404 WRITE (numout,*) ' iper = 13' 405 ENDIF 406 iperm1 = 12 407 iper =13 408 ENDIF 409 ! 410 ! we have to prepare a NEW READ of DATA 411 ! 412 ! swap from record 2 to 1 413 ! 414 udta(:,:,:,1) = udta(:,:,:,2) 415 vdta(:,:,:,1) = vdta(:,:,:,2) 416 wdta(:,:,:,1) = wdta(:,:,:,2) 417 #if defined key_trc_diatrd 418 hdivdta(:,:,:,1)=hdivdta(:,:,:,2) 419 #endif 420 avtdta(:,:,:,1) = avtdta(:,:,:,2) 421 tdta(:,:,:,1) = tdta(:,:,:,2) 422 sdta(:,:,:,1) = sdta(:,:,:,2) 271 IF( iperm1 /= ndyn1 ) THEN 272 273 IF( iperm1 == 0. ) THEN 274 IF (lwp) THEN 275 WRITE (numout,*) ' dynamic file is not periodic with periodic interpolation' 276 WRITE (numout,*) ' we take the last value for the last period ' 277 WRITE (numout,*) ' iperm1 = 12, iper = 13 ' 278 ENDIF 279 iperm1 = 12 280 iper = 13 281 ENDIF 282 ! 283 ! We have to prepare a new read of data : swap from record 2 to 1 284 ! 285 CALL swap_dyn_data 286 287 iswap = 1 ! indicates swap 288 289 CALL dynrea( kt, iper ) ! data read for the iper period 290 423 291 #if defined key_ldfslp 424 uslpdta(:,:,:,1) = uslpdta(:,:,:,2) 425 vslpdta(:,:,:,1) = vslpdta(:,:,:,2) 426 wslpidta(:,:,:,1) = wslpidta(:,:,:,2) 427 wslpjdta(:,:,:,1) = wslpjdta(:,:,:,2) 428 #endif 429 flxdta(:,:,:,1) = flxdta(:,:,:,2) 430 zmxldta(:,:,1) = zmxldta(:,:,2) 431 432 #if ! defined key_off_degrad 433 434 # if defined key_traldf_c2d 435 ahtwdta(:,:,1)= ahtwdta(:,:,2) 436 # if defined key_trcldf_eiv 437 aeiwdta(:,:,1)= aeiwdta(:,:,2) 438 # endif 439 # endif 440 441 #else 442 ahtudta(:,:,:,1) = ahtudta(:,:,:,2) 443 ahtvdta(:,:,:,1) = ahtvdta(:,:,:,2) 444 ahtwdta(:,:,:,1) = ahtwdta(:,:,:,2) 445 # if defined key_trcldf_eiv 446 aeiudta(:,:,:,1) = aeiudta(:,:,:,2) 447 aeivdta(:,:,:,1) = aeivdta(:,:,:,2) 448 aeiwdta(:,:,:,1) = aeiwdta(:,:,:,2) 449 # endif 450 451 #endif 452 453 #if defined key_trcbbl_dif || defined key_trcbbl_adv 454 bblxdta(:,:,1) = bblxdta(:,:,2) 455 bblydta(:,:,1) = bblydta(:,:,2) 456 #endif 457 ! 458 ! indicates a swap 459 ! 460 iswap = 1 461 ! 462 ! READ DATA for the iper period 463 ! 464 CALL dynrea( kt, iper ) 465 ! 466 ! Computes wdta (and slopes if key_trahdfiso) 467 ! 468 tn(:,:,:)=tdta(:,:,:,2) 469 sn(:,:,:)=sdta(:,:,:,2) 470 avt(:,:,:)=avtdta(:,:,:,2) 471 472 #if defined key_ldfslp 473 CALL eos( tn, sn, rhd, rhop ) ! Time-filtered in situ density 474 CALL bn2( tn, sn, rn2 ) ! before Brunt-Vaisala frequency 475 IF( ln_zps ) & 476 CALL zps_hde( kt, tn , sn , rhd, & ! Partial steps: before Horizontal DErivative 477 gtu, gsu, gru, & ! of t, s, rd at the bottom ocean level 478 gtv, gsv, grv ) 479 CALL zdf_mxl( kt ) ! mixed layer depth 480 CALL ldf_slp( kt, rhd, rn2 ) 481 482 uslpdta(:,:,:,2)=uslp(:,:,:) 483 vslpdta(:,:,:,2)=vslp(:,:,:) 484 wslpidta(:,:,:,2)=wslpi(:,:,:) 485 wslpjdta(:,:,:,2)=wslpj(:,:,:) 486 #endif 487 ! 488 ! store the information of the period read 489 ! 490 ndyn1 = ndyn2 491 ndyn2 = iper 492 ! 493 ! we have READ another period of DATA ! 494 IF (lwp) THEN 495 WRITE (numout,*) ' dynamics DATA READ for the period ndyn1 =',ndyn1 496 WRITE (numout,*) ' and the period ndyn2 = ',ndyn2 497 WRITE (numout,*) ' time step is :',kt 498 END IF 499 500 END IF 501 502 ! 503 ! compute the DATA at the given time step 504 ! 505 IF ( nsptint == 0 ) THEN 506 ! 507 ! no spatial interpolation 508 ! 509 ! DATA are probably correct 510 ! we have to initialize DATA IF we have changed the period 511 ! 512 IF (iswap.eq.1) THEN 513 ! 514 ! initialize now fields with the NEW DATA READ 515 ! 516 un(:,:,:)=udta(:,:,:,2) 517 vn(:,:,:)=vdta(:,:,:,2) 518 wn(:,:,:)=wdta(:,:,:,2) 519 hdivn(:,:,:)=hdivdta(:,:,:,2) 520 #if defined key_trc_zdfddm 521 avs(:,:,:)=avtdta(:,:,:,2) 522 #endif 523 avt(:,:,:)=avtdta(:,:,:,2) 524 tn(:,:,:)=tdta(:,:,:,2) 525 sn(:,:,:)=sdta(:,:,:,2) 526 #if defined key_ldfslp 527 uslp(:,:,:)=uslpdta(:,:,:,2) 528 vslp(:,:,:)=vslpdta(:,:,:,2) 529 wslpi(:,:,:)=wslpidta(:,:,:,2) 530 wslpj(:,:,:)=wslpjdta(:,:,:,2) 531 #endif 532 flx(:,:,:) = flxdta(:,:,:,2) 533 hmld(:,:)=zmxldta(:,:,2) 534 #if ! defined key_off_degrad 535 536 # if defined key_traldf_c2d 537 ahtwdta(:,:,1)= ahtwdta(:,:,2) 538 # if defined key_trcldf_eiv 539 aeiwdta(:,:,1)= aeiwdta(:,:,2) 540 # endif 541 # endif 542 543 #else 544 ahtudta(:,:,:,1) = ahtudta(:,:,:,2) 545 ahtvdta(:,:,:,1) = ahtvdta(:,:,:,2) 546 ahtwdta(:,:,:,1) = ahtwdta(:,:,:,2) 547 # if defined key_trcldf_eiv 548 aeiudta(:,:,:,1) = aeiudta(:,:,:,2) 549 aeivdta(:,:,:,1) = aeivdta(:,:,:,2) 550 aeiwdta(:,:,:,1) = aeiwdta(:,:,:,2) 551 # endif 552 553 #endif 554 555 #if defined key_trcbbl_dif || defined key_trcbbl_adv 556 bblx(:,:)=bblxdta(:,:,2) 557 bbly(:,:)=bblydta(:,:,2) 558 #endif 559 ! 560 ! keep needed fluxes 561 ! 562 wndm(:,:) = flx(:,:,jpwind) 563 fr_i(:,:) = flx(:,:,jpice) 564 emp(:,:) = flx(:,:,jpemp) 565 emps(:,:) = emp(:,:) 566 qsr(:,:) = flx(:,:,jpqsr) 567 568 END IF 569 570 ELSE 571 IF ( nsptint == 1 ) THEN 572 ! 573 ! linear interpolation 574 ! 575 ! initialize now fields with a linear interpolation 576 ! 577 un(:,:,:) = zweighm1 * udta(:,:,:,1) + zweigh * udta(:,:,:,2) 578 vn(:,:,:) = zweighm1 * vdta(:,:,:,1) + zweigh * vdta(:,:,:,2) 579 wn(:,:,:) = zweighm1 * wdta(:,:,:,1) + zweigh * wdta(:,:,:,2) 580 hdivn(:,:,:) = zweighm1 * hdivdta(:,:,:,1) + zweigh * hdivdta(:,:,:,2) 581 #if defined key_zdfddm 582 avs(:,:,:)= zweighm1 * avtdta(:,:,:,1) + zweigh * avtdta(:,:,:,2) 583 #endif 584 avt(:,:,:)= zweighm1 * avtdta(:,:,:,1) + zweigh * avtdta(:,:,:,2) 585 tn(:,:,:) = zweighm1 * tdta(:,:,:,1) + zweigh * tdta(:,:,:,2) 586 sn(:,:,:) = zweighm1 * sdta(:,:,:,1) + zweigh * sdta(:,:,:,2) 587 588 589 #if defined key_ldfslp 590 uslp(:,:,:) = zweighm1 * uslpdta(:,:,:,1) + zweigh * uslpdta(:,:,:,2) 591 vslp(:,:,:) = zweighm1 * vslpdta(:,:,:,1) + zweigh * vslpdta(:,:,:,2) 592 wslpi(:,:,:) = zweighm1 * wslpidta(:,:,:,1) + zweigh * wslpidta(:,:,:,2) 593 wslpj(:,:,:) = zweighm1 * wslpjdta(:,:,:,1) + zweigh * wslpjdta(:,:,:,2) 594 #endif 595 flx(:,:,:) = zweighm1 * flxdta(:,:,:,1) + zweigh * flxdta(:,:,:,2) 596 hmld(:,:) = zweighm1 * zmxldta(:,:,1) + zweigh * zmxldta(:,:,2) 597 #if ! defined key_off_degrad 598 599 # if defined key_traldf_c2d 600 ahtw(:,:) = zweighm1 * ahtwdta(:,:,1) + zweigh * ahtwdta(:,:,2) 601 # if defined key_trcldf_eiv 602 aeiw(:,:) = zweighm1 * aeiwdta(:,:,1) + zweigh * aeiwdta(:,:,2) 603 # endif 604 # endif 605 606 #else 607 ahtu(:,:,:) = zweighm1 * ahtudta(:,:,:,1) + zweigh * ahtudta(:,:,:,2) 608 ahtv(:,:,:) = zweighm1 * ahtvdta(:,:,:,1) + zweigh * ahtvdta(:,:,:,2) 609 ahtw(:,:,:) = zweighm1 * ahtwdta(:,:,:,1) + zweigh * ahtwdta(:,:,:,2) 610 # if defined key_trcldf_eiv 611 aeiu(:,:,:) = zweighm1 * aeiudta(:,:,:,1) + zweigh * aeiudta(:,:,:,2) 612 aeiv(:,:,:) = zweighm1 * aeivdta(:,:,:,1) + zweigh * aeivdta(:,:,:,2) 613 aeiw(:,:,:) = zweighm1 * aeiwdta(:,:,:,1) + zweigh * aeiwdta(:,:,:,2) 614 # endif 615 616 #endif 617 618 #if defined key_trcbbl_dif || defined key_trcbbl_adv 619 bblx(:,:) = zweighm1 * bblxdta(:,:,1) + zweigh * bblxdta(:,:,2) 620 bbly(:,:) = zweighm1 * bblydta(:,:,1) + zweigh * bblydta(:,:,2) 621 #endif 622 ! 623 ! keep needed fluxes 624 ! 625 wndm(:,:) = flx(:,:,jpwind) 626 fr_i(:,:) = flx(:,:,jpice) 627 emp(:,:) = flx(:,:,jpemp) 628 emps(:,:) = emp(:,:) 629 qsr(:,:) = flx(:,:,jpqsr) 630 ! 631 ! other interpolation 632 ! 633 ELSE 634 635 WRITE (numout,*) ' this kind of interpolation don t EXIST' 636 WRITE (numout,*) ' at the moment. we STOP ' 637 STOP 'dtadyn' 638 639 END IF 640 292 ! Computes slopes 293 ! Caution : here tn, sn and avt are used as workspace 294 tn (:,:,:) = tdta (:,:,:,2) 295 sn (:,:,:) = sdta (:,:,:,2) 296 avt(:,:,:) = avtdta(:,:,:,2) 297 298 CALL eos( tn, sn, rhd, rhop ) ! Time-filtered in situ density 299 CALL bn2( tn, sn, rn2 ) ! before Brunt-Vaisala frequency 300 IF( ln_zps ) & 301 & CALL zps_hde( kt, tn , sn , rhd, & ! Partial steps: before Horizontal DErivative 302 & gtu, gsu, gru, & ! of t, s, rd at the bottom ocean level 303 & gtv, gsv, grv ) 304 CALL zdf_mxl( kt ) ! mixed layer depth 305 CALL ldf_slp( kt, rhd, rn2 ) 306 307 uslpdta (:,:,:,2) = uslp (:,:,:) 308 vslpdta (:,:,:,2) = vslp (:,:,:) 309 wslpidta(:,:,:,2) = wslpi(:,:,:) 310 wslpjdta(:,:,:,2) = wslpj(:,:,:) 311 #endif 312 313 ! store the information of the period read 314 ndyn1 = ndyn2 315 ndyn2 = iper 316 317 IF (lwp) THEN 318 WRITE (numout,*) ' dynamics data read for the period ndyn1 =',ndyn1, & 319 & ' and for the period ndyn2 = ',ndyn2 320 WRITE (numout,*) ' time step is : ', kt 321 END IF 322 ! 641 323 END IF 642 324 ! 643 ! lb in any case, we need rhop 644 ! 325 ! Compute the data at the given time step 326 !---------------------------------------- 327 328 IF( nsptint == 0 ) THEN 329 ! No spatial interpolation, data are probably correct 330 ! We have to initialize data if we have changed the period 331 CALL assign_dyn_data 332 ELSE IF( nsptint == 1 ) THEN 333 ! linear interpolation 334 CALL linear_interp_dyn_data( zweigh ) 335 ELSE 336 ! other interpolation 337 WRITE (numout,*) ' this kind of interpolation do not exist at the moment : we stop' 338 STOP 'dtadyn' 339 END IF 340 341 ! In any case, we need rhop 645 342 CALL eos( tn, sn, rhd, rhop ) 646 343 647 344 #if ! defined key_off_degrad && defined key_traldf_c2d 648 345 ! In case of 2D varying coefficients, we need aeiv and aeiu 649 346 IF( lk_traldf_eiv ) CALL ldf_eiv( kt ) ! eddy induced velocity coefficient 650 347 #endif 651 348 652 349 END SUBROUTINE dta_dyn 653 350 … … 672 369 INTEGER, INTENT( in ) :: kt, kenr ! time index 673 370 !! * Local declarations 674 INTEGER :: j i, jj, jk, jkenr371 INTEGER :: jkenr 675 372 676 373 REAL(wp), DIMENSION(jpi,jpj,jpk) :: & … … 679 376 680 377 REAL(wp), DIMENSION(jpi,jpj) :: & 681 zemp, zqsr, zmld, zice, zwspd 378 zemp, zqsr, zmld, zice, zwspd, & 379 ztaux, ztauy 682 380 #if defined key_trcbbl_dif || defined key_trcbbl_adv 683 381 REAL(wp), DIMENSION(jpi,jpj) :: zbblx, zbbly 684 382 #endif 685 383 686 #if ! defined key_off_degrad 687 # if defined key_traldf_c2d 384 #if ! defined key_off_degrad && defined key_traldf_c2d 688 385 REAL(wp), DIMENSION(jpi,jpj) :: zahtw 689 386 # if defined key_trcldf_eiv 690 387 REAL(wp), DIMENSION(jpi,jpj) :: zaeiw 691 # 692 # 693 694 # else388 # endif 389 #endif 390 391 #if defined key_off_degrad 695 392 REAL(wp), DIMENSION(jpi,jpj,jpk) :: & 696 393 zahtu, zahtv, zahtw ! Lateral diffusivity … … 699 396 zaeiu, zaeiv, zaeiw ! G&M coefficient 700 397 # endif 701 702 #endif 703 704 # if defined key_diaeiv 705 !! GM Velocity : to be used latter 706 REAL(wp), DIMENSION(jpi,jpj,jpk) :: & 707 zeivu, zeivv, zeivw 708 # endif 709 710 CHARACTER(len=45) :: & 711 clname_t = 'dyna_grid_T.nc', & 712 clname_u = 'dyna_grid_U.nc', & 713 clname_v = 'dyna_grid_V.nc', & 714 clname_w = 'dyna_grid_W.nc' 398 #endif 399 715 400 !--------------------------------------------------------------- 716 401 ! 0. Initialization … … 733 418 734 419 IF( kt == nit000 .AND. nlecoff == 0 ) THEN 735 736 420 nlecoff = 1 737 738 CALL iom_open ( clname_t, numfl_t ) 739 CALL iom_open ( clname_u, numfl_u ) 740 CALL iom_open ( clname_v, numfl_v ) 741 CALL iom_open ( clname_w, numfl_w ) 742 421 CALL iom_open ( cfile_grid_T, numfl_t ) 422 CALL iom_open ( cfile_grid_U, numfl_u ) 423 CALL iom_open ( cfile_grid_V, numfl_v ) 424 CALL iom_open ( cfile_grid_W, numfl_w ) 743 425 ENDIF 744 426 745 427 ! file grid-T 746 428 !--------------- 747 CALL iom_get ( numfl_t, jpdom_data, 'votemper', zt (:,:,:), jkenr ) 748 CALL iom_get ( numfl_t, jpdom_data, 'vosaline', zs (:,:,:), jkenr ) 749 CALL iom_get ( numfl_t, jpdom_data, 'somixhgt', zmld (:,: ), jkenr ) 750 CALL iom_get ( numfl_t, jpdom_data, 'sowaflcd', zemp (:,: ), jkenr ) 751 CALL iom_get ( numfl_t, jpdom_data, 'soshfldo', zqsr (:,: ), jkenr ) 752 CALL iom_get ( numfl_t, jpdom_data, 'soicecov', zice (:,: ), jkenr ) 753 CALL iom_get ( numfl_t, jpdom_data, 'sowindsp', zwspd(:,: ), jkenr ) 754 755 ! file grid-U 756 !--------------- 757 CALL iom_get ( numfl_u, jpdom_data, 'vozocrtx', zu (:,:,:), jkenr ) 758 #if defined key_trcbbl_dif || defined key_trcbbl_adv 759 CALL iom_get ( numfl_u, jpdom_data, 'sobblcox', zbblx(:,: ), jkenr ) 760 #endif 761 762 #if defined key_diaeiv 763 !! GM Velocity : to be used latter 764 CALL iom_get ( numfl_u, jpdom_data, 'vozoeivu', zeivu(:,:,:), jkenr ) 765 #endif 766 767 # if defined key_off_degrad 768 CALL iom_get ( numfl_u, jpdom_data, 'vozoahtu', zahtu(:,:,:), jkenr ) 769 # if defined key_trcldf_eiv 770 CALL iom_get ( numfl_u, jpdom_data, 'vozoaeiu', zaeiu(:,:,:), jkenr ) 771 # endif 772 #endif 773 774 ! file grid-V 775 !--------------- 776 CALL iom_get ( numfl_v, jpdom_data, 'vomecrty', zv (:,:,:), jkenr ) 777 #if defined key_trcbbl_dif || defined key_trcbbl_adv 778 CALL iom_get ( numfl_v, jpdom_data, 'sobblcoy', zbbly(:,: ), jkenr ) 779 #endif 780 781 #if defined key_diaeiv 782 !! GM Velocity : to be used latter 783 CALL iom_get ( numfl_v, jpdom_data, 'vomeeivv', zeivv(:,:,:), jkenr ) 429 CALL iom_get( numfl_t, jpdom_data, 'votemper', zt (:,:,:), jkenr ) 430 CALL iom_get( numfl_t, jpdom_data, 'vosaline', zs (:,:,:), jkenr ) 431 CALL iom_get( numfl_t, jpdom_data, 'somixhgt', zmld (:,: ), jkenr ) 432 CALL iom_get( numfl_t, jpdom_data, 'sowaflcd', zemp (:,: ), jkenr ) 433 CALL iom_get( numfl_t, jpdom_data, 'soshfldo', zqsr (:,: ), jkenr ) 434 CALL iom_get( numfl_t, jpdom_data, 'soicecov', zice (:,: ), jkenr ) 435 IF( iom_varid( numfl_t, 'sowindsp', ldstop = .FALSE. ) > 0 ) THEN 436 CALL iom_get( numfl_t, jpdom_data, 'sowindsp', zwspd(:,:), jkenr ) 437 ELSE 438 CALL iom_get( numfl_u, jpdom_data, 'sozotaux', ztaux(:,:), jkenr ) 439 CALL iom_get( numfl_v, jpdom_data, 'sometauy', ztauy(:,:), jkenr ) 440 CALL tau2wnd( ztaux, ztauy, zwspd ) 441 ENDIF 442 ! files grid-U / grid_V 443 CALL iom_get( numfl_u, jpdom_data, 'vozocrtx', zu (:,:,:), jkenr ) 444 CALL iom_get( numfl_v, jpdom_data, 'vomecrty', zv (:,:,:), jkenr ) 445 446 #if defined key_trcbbl_dif || defined key_trcbbl_adv 447 IF( iom_varid( numfl_u, 'sobblcox', ldstop = .FALSE. ) > 0 .AND. & 448 & iom_varid( numfl_v, 'sobblcoy', ldstop = .FALSE. ) > 0 ) THEN 449 CALL iom_get( numfl_u, jpdom_data, 'sobblcox', zbblx(:,:), jkenr ) 450 CALL iom_get( numfl_v, jpdom_data, 'sobblcoy', zbbly(:,:), jkenr ) 451 ELSE 452 CALL bbl_sign( zt, zs, zbblx, zbbly ) 453 ENDIF 454 #endif 455 456 ! file grid-W 457 !! CALL iom_get ( numfl_w, jpdom_data, 'vovecrtz', zw (:,:,:), jkenr ) 458 ! Computation of vertical velocity using horizontal divergence 459 CALL wzv( zu, zv, zw, zhdiv ) 460 461 # if defined key_zdfddm 462 CALL iom_get( numfl_w, jpdom_data, 'voddmavs', zavt (:,:,:), jkenr ) 463 #else 464 CALL iom_get( numfl_w, jpdom_data, 'votkeavt', zavt (:,:,:), jkenr ) 465 #endif 466 467 #if ! defined key_off_degrad && defined key_traldf_c2d 468 CALL iom_get( numfl_w, jpdom_data, 'soleahtw', zahtw (:,: ), jkenr ) 469 # if defined key_trcldf_eiv 470 CALL iom_get( numfl_w, jpdom_data, 'soleaeiw', zaeiw (:,: ), jkenr ) 471 # endif 784 472 #endif 785 473 786 474 #if defined key_off_degrad 787 CALL iom_get ( numfl_v, jpdom_data, 'vomeahtv', zahtv(:,:,:), jkenr ) 788 # if defined key_trcldf_eiv 789 CALL iom_get ( numfl_v, jpdom_data, 'vomeaeiv', zaeiv(:,:,:), jkenr ) 790 # endif 791 #endif 792 793 ! file grid-W 794 !--------------- 795 !! CALL iom_get ( numfl_w, jpdom_data, 'vovecrtz', zw (:,:,:), jkenr ) 796 # if defined key_zdfddm 797 CALL iom_get ( numfl_w, jpdom_data, 'voddmavs', zavt (:,:,:), jkenr ) 798 #else 799 CALL iom_get ( numfl_w, jpdom_data, 'votkeavt', zavt (:,:,:), jkenr ) 800 #endif 801 802 # if defined key_diaeiv 803 !! GM Velocity : to be used latter 804 CALL iom_get ( numfl_w, jpdom_data, 'voveeivw', zeivw(:,:,:), jkenr ) 805 #endif 806 807 #if ! defined key_off_degrad 808 # if defined key_traldf_c2d 809 CALL iom_get ( numfl_w, jpdom_data, 'soleahtw', zahtw (:,: ), jkenr ) 810 # if defined key_trcldf_eiv 811 CALL iom_get ( numfl_w, jpdom_data, 'soleaeiw', zaeiw (:,: ), jkenr ) 812 # endif 813 # endif 814 #else 815 !! degradation-integration 816 CALL iom_get ( numfl_w, jpdom_data, 'voveahtw', zahtw(:,:,:), jkenr ) 817 # if defined key_trcldf_eiv 818 CALL iom_get ( numfl_w, jpdom_data, 'voveaeiw', zaeiw(:,:,:), jkenr ) 819 # endif 475 CALL iom_get( numfl_u, jpdom_data, 'vozoahtu', zahtu(:,:,:), jkenr ) 476 CALL iom_get( numfl_v, jpdom_data, 'vomeahtv', zahtv(:,:,:), jkenr ) 477 CALL iom_get( numfl_w, jpdom_data, 'voveahtw', zahtw(:,:,:), jkenr ) 478 # if defined key_trcldf_eiv 479 CALL iom_get( numfl_u, jpdom_data, 'vozoaeiu', zaeiu(:,:,:), jkenr ) 480 CALL iom_get( numfl_v, jpdom_data, 'vomeaeiv', zaeiv(:,:,:), jkenr ) 481 CALL iom_get( numfl_w, jpdom_data, 'voveaeiw', zaeiw(:,:,:), jkenr ) 482 # endif 820 483 #endif 821 484 822 485 udta(:,:,:,2) = zu(:,:,:) * umask(:,:,:) 823 vdta(:,:,:,2) = zv(:,:,:) * vmask(:,:,:) 824 !! wdta(:,:,:,2) = zw(:,:,:) * tmask(:,:,:) 825 826 827 ! Computation of vertical velocity using horizontal divergence 828 zhdiv(:,:,:) = 0. 829 DO jk = 1, jpkm1 830 DO jj = 2, jpjm1 831 DO ji = fs_2, fs_jpim1 ! vector opt. 832 #if defined key_zco 833 zhdiv(ji,jj,jk) = ( e2u(ji,jj) * udta(ji,jj,jk,2) - e2u(ji-1,jj ) * udta(ji-1,jj ,jk,2) & 834 & + e1v(ji,jj) * vdta(ji,jj,jk,2) - e1v(ji ,jj-1) * vdta(ji ,jj-1,jk,2) ) & 835 & / ( e1t(ji,jj) * e2t(ji,jj) ) 836 #else 837 zhdiv(ji,jj,jk) = & 838 ( 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) & 839 + 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) ) & 840 / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) 841 #endif 842 END DO 843 END DO 844 ENDDO 845 846 ! Lateral boundary conditions on hdiv 847 CALL lbc_lnk( zhdiv, 'T', 1. ) 848 849 850 ! computation of vertical velocity from the bottom 851 zw(:,:,jpk) = 0. 852 DO jk = jpkm1, 1, -1 853 zw(:,:,jk) = zw(:,:,jk+1) - fse3t(:,:,jk) * zhdiv(:,:,jk) 854 END DO 855 486 vdta(:,:,:,2) = zv(:,:,:) * vmask(:,:,:) 856 487 wdta(:,:,:,2) = zw(:,:,:) * tmask(:,:,:) 488 489 #if defined key_trc_diatrd 857 490 hdivdta(:,:,:,2) = zhdiv(:,:,:) * tmask(:,:,:) 858 859 tdta(:,:,:,2) = zt(:,:,:) * tmask(:,:,:) 860 sdta(:,:,:,2) = zs(:,:,:) * tmask(:,:,:) 491 #endif 492 493 tdta(:,:,:,2) = zt (:,:,:) * tmask(:,:,:) 494 sdta(:,:,:,2) = zs (:,:,:) * tmask(:,:,:) 861 495 avtdta(:,:,:,2) = zavt(:,:,:) * tmask(:,:,:) 496 862 497 #if ! defined key_off_degrad && defined key_traldf_c2d 863 498 ahtwdta(:,:,2) = zahtw(:,:) * tmask(:,:,1) … … 878 513 #endif 879 514 515 ! fluxes 880 516 ! 881 ! flux : 882 ! 883 flxdta(:,:,jpwind,2) = zwspd(:,:) * tmask(:,:,1) 884 flxdta(:,:,jpice,2) = MIN( 1., zice(:,:) ) * tmask(:,:,1) 885 flxdta(:,:,jpemp,2) = zemp(:,:) * tmask(:,:,1) 886 flxdta(:,:,jpqsr,2) = zqsr(:,:) * tmask(:,:,1) 887 zmxldta(:,:,2) = zmld(:,:) * tmask(:,:,1) 517 wspddta(:,:,2) = zwspd(:,:) * tmask(:,:,1) 518 frlddta(:,:,2) = MIN( 1., zice(:,:) ) * tmask(:,:,1) 519 empdta (:,:,2) = zemp(:,:) * tmask(:,:,1) 520 qsrdta (:,:,2) = zqsr(:,:) * tmask(:,:,1) 521 hmlddta(:,:,2) = zmld(:,:) * tmask(:,:,1) 888 522 889 523 #if defined key_trcbbl_dif || defined key_trcbbl_adv … … 893 527 WHERE( bblxdta(:,:,2) > 2. ) bblxdta(:,:,2) = 0. 894 528 WHERE( bblydta(:,:,2) > 2. ) bblydta(:,:,2) = 0. 895 896 529 #endif 897 530 … … 905 538 END SUBROUTINE dynrea 906 539 907 540 SUBROUTINE dta_dyn_init 541 !!---------------------------------------------------------------------- 542 !! *** ROUTINE dta_dyn_init *** 543 !! 544 !! ** Purpose : initializations of parameters for the interpolation 545 !! 546 !! ** Method : 547 !! 548 !! History : 549 !! ! original : 92-01 (M. Imbard: sub domain) 550 !! ! 98-04 (L.Bopp MA Foujols: slopes for isopyc.) 551 !! ! 98-05 (L. Bopp read output of coupled run) 552 !! ! 05-03 (O. Aumont and A. El Moussaoui) F90 553 !!---------------------------------------------------------------------- 554 !! * Modules used 555 556 !! * Local declarations 557 558 559 NAMELIST/namdyn/ ndtadyn, ndtatot, nsptint, nficdyn, lperdyn, & 560 & cfile_grid_T, cfile_grid_U, cfile_grid_V, cfile_grid_W 561 !!---------------------------------------------------------------------- 562 563 ! Define the dynamical input parameters 564 ! ====================================== 565 566 ! Read Namelist namdyn : Lateral physics on tracers 567 REWIND( numnam ) 568 READ ( numnam, namdyn ) 569 570 IF(lwp) THEN 571 WRITE(numout,*) 572 WRITE(numout,*) 'namdyn : offline dynamical selection' 573 WRITE(numout,*) '~~~~~~~' 574 WRITE(numout,*) ' Namelist namdyn : set parameters for the lecture of the dynamical fields' 575 WRITE(numout,*) 576 WRITE(numout,*) ' number of elements in the FILE for a year ndtadyn = ' , ndtadyn 577 WRITE(numout,*) ' total number of elements in the FILE ndtatot = ' , ndtatot 578 WRITE(numout,*) ' type of interpolation nsptint = ' , nsptint 579 WRITE(numout,*) ' number of dynamics FILE nficdyn = ' , nficdyn 580 WRITE(numout,*) ' loop on the same FILE lperdyn = ' , lperdyn 581 WRITE(numout,*) ' ' 582 WRITE(numout,*) ' name of grid_T file cfile_grid_T = ', TRIM(cfile_grid_T) 583 WRITE(numout,*) ' name of grid_U file cfile_grid_U = ', TRIM(cfile_grid_U) 584 WRITE(numout,*) ' name of grid_V file cfile_grid_V = ', TRIM(cfile_grid_V) 585 WRITE(numout,*) ' name of grid_W file cfile_grid_W = ', TRIM(cfile_grid_W) 586 WRITE(numout,*) ' ' 587 ENDIF 588 589 END SUBROUTINE dta_dyn_init 590 591 SUBROUTINE wzv( pu, pv, pw, phdiv ) 592 !!---------------------------------------------------------------------- 593 !! *** ROUTINE wzv *** 594 !! 595 !! ** Purpose : Compute the now vertical velocity after the array swap 596 !! 597 !! ** Method : 598 !! ** Method : - Divergence: 599 !! - compute the now divergence given by : 600 !! * z-coordinate 601 !! hdiv = 1/(e1t*e2t) [ di(e2u u) + dj(e1v v) ] 602 !! - Using the incompressibility hypothesis, the vertical 603 !! velocity is computed by integrating the horizontal divergence 604 !! from the bottom to the surface. 605 !! The boundary conditions are w=0 at the bottom (no flux) and, 606 !! in regid-lid case, w=0 at the sea surface. 607 !! 608 !! 609 !! History : 610 !! 9.0 ! 02-07 (G. Madec) Vector optimization 611 !!---------------------------------------------------------------------- 612 !! * Arguments 613 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( in ) :: pu, pv !: horizontal velocities 614 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( out ) :: pw !: verticla velocity 615 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( inout ) :: phdiv !: horizontal divergence 616 617 !! * Local declarations 618 INTEGER :: ji, jj, jk 619 REAL(wp) :: zu, zu1, zv, zv1, zet 620 621 622 ! Computation of vertical velocity using horizontal divergence 623 phdiv(:,:,:) = 0. 624 DO jk = 1, jpkm1 625 DO jj = 2, jpjm1 626 DO ji = fs_2, fs_jpim1 ! vector opt. 627 #if defined key_zco 628 zu = pu(ji ,jj ,jk) * umask(ji ,jj ,jk) * e2u(ji ,jj ) 629 zu1 = pu(ji-1,jj ,jk) * umask(ji-1,jj ,jk) * e2u(ji-1,jj ) 630 zv = pv(ji ,jj ,jk) * vmask(ji ,jj ,jk) * e1v(ji ,jj ) 631 zv1 = pv(ji ,jj-1,jk) * vmask(ji ,jj-1,jk) * e1v(ji ,jj-1) 632 zet = 1. / ( e1t(ji,jj) * e2t(ji,jj) ) 633 #else 634 zu = pu(ji ,jj ,jk) * umask(ji ,jj ,jk) * e2u(ji ,jj ) * fse3u(ji ,jj ,jk) 635 zu1 = pu(ji-1,jj ,jk) * umask(ji-1,jj ,jk) * e2u(ji-1,jj ) * fse3u(ji-1,jj ,jk) 636 zv = pv(ji ,jj ,jk) * vmask(ji ,jj ,jk) * e1v(ji ,jj ) * fse3v(ji ,jj ,jk) 637 zv1 = pv(ji ,jj-1,jk) * vmask(ji ,jj-1,jk) * e1v(ji ,jj-1) * fse3v(ji ,jj-1,jk) 638 zet = 1. / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) 639 #endif 640 phdiv(ji,jj,jk) = ( zu - zu1 + zv - zv1 ) * zet 641 END DO 642 END DO 643 ENDDO 644 645 ! Lateral boundary conditions on phdiv 646 CALL lbc_lnk( phdiv, 'T', 1. ) 647 648 649 ! computation of vertical velocity from the bottom 650 pw(:,:,jpk) = 0. 651 DO jk = jpkm1, 1, -1 652 pw(:,:,jk) = pw(:,:,jk+1) - fse3t(:,:,jk) * phdiv(:,:,jk) 653 END DO 654 655 END SUBROUTINE wzv 656 657 SUBROUTINE tau2wnd( ptaux, ptauy, pwspd ) 658 !!--------------------------------------------------------------------- 659 !! *** ROUTINE sbc_tau2wnd *** 660 !! 661 !! ** Purpose : Estimation of wind speed as a function of wind stress 662 !! 663 !! ** Method : |tau|=rhoa*Cd*|U|^2 664 !!--------------------------------------------------------------------- 665 !! * Arguments 666 REAL(wp), DIMENSION(jpi,jpj), INTENT( in ) :: & 667 ptaux, ptauy !: wind stress in i-j direction resp. 668 REAL(wp), DIMENSION(jpi,jpj), INTENT( out ) :: & 669 pwspd !: wind speed 670 REAL(wp) :: zrhoa = 1.22 ! Air density kg/m3 671 REAL(wp) :: zcdrag = 1.5e-3 ! drag coefficient 672 REAL(wp) :: ztx, zty, ztau, zcoef ! temporary variables 673 INTEGER :: ji, jj ! dummy indices 674 !!--------------------------------------------------------------------- 675 zcoef = 0.5 / ( zrhoa * zcdrag ) 676 !CDIR NOVERRCHK 677 DO jj = 2, jpjm1 678 !CDIR NOVERRCHK 679 DO ji = fs_2, fs_jpim1 ! vect. opt. 680 ztx = ptaux(ji-1,jj ) + ptaux(ji,jj) 681 zty = ptauy(ji ,jj-1) + ptauy(ji,jj) 682 ztau = SQRT( ztx * ztx + zty * zty ) 683 pwspd(ji,jj) = SQRT ( ztau * zcoef ) * tmask(ji,jj,1) 684 END DO 685 END DO 686 CALL lbc_lnk( pwspd(:,:) , 'T', 1. ) 687 688 END SUBROUTINE tau2wnd 689 690 #if defined key_trcbbl_dif || defined key_trcbbl_adv 691 692 SUBROUTINE bbl_sign( ptn, psn, pbblx, pbbly ) 693 !!---------------------------------------------------------------------- 694 !! *** ROUTINE bbl_sign *** 695 !! 696 !! ** Purpose : Compute the sign of local gradient of density multiplied by the slope 697 !! along the bottom slope gradient : grad( rho) * grad(h) 698 !! Need to compute the diffusive bottom boundary layer 699 !! 700 !! ** Method : When the product grad( rho) * grad(h) < 0 (where grad 701 !! is an along bottom slope gradient) an additional lateral diffu- 702 !! sive trend along the bottom slope is added to the general tracer 703 !! trend, otherwise nothing is done. See trcbbl.F90 704 !! 705 !! 706 !! History : 707 !! 9.0 ! 02-07 (G. Madec) Vector optimization 708 !!---------------------------------------------------------------------- 709 !! * Arguments 710 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( in ) :: & 711 ptn , & !: temperature 712 psn !: salinity 713 REAL(wp), DIMENSION(jpi,jpj), INTENT( out ) :: & 714 pbblx , pbbly !: sign of bbl in i-j direction resp. 715 716 !! * Local declarations 717 INTEGER :: ji, jj ! dummy loop indices 718 INTEGER :: ik 719 REAL(wp) :: & 720 ztx, zsx, zhx, zalbetx, zgdrhox, & ! temporary scalars 721 zty, zsy, zhy, zalbety, zgdrhoy 722 REAL(wp), DIMENSION(jpi,jpj) :: & 723 ztnb, zsnb, zdep 724 REAL(wp) :: fsalbt, pft, pfs, pfh ! statement function 725 !!---------------------------------------------------------------------- 726 ! ratio alpha/beta 727 ! ================ 728 ! fsalbt: ratio of thermal over saline expension coefficients 729 ! pft : potential temperature in degrees celcius 730 ! pfs : salinity anomaly (s-35) in psu 731 ! pfh : depth in meters 732 733 fsalbt( pft, pfs, pfh ) = & 734 ( ( ( -0.255019e-07 * pft + 0.298357e-05 ) * pft & 735 - 0.203814e-03 ) * pft & 736 + 0.170907e-01 ) * pft & 737 + 0.665157e-01 & 738 +(-0.678662e-05 * pfs - 0.846960e-04 * pft + 0.378110e-02 ) * pfs & 739 + ( ( - 0.302285e-13 * pfh & 740 - 0.251520e-11 * pfs & 741 + 0.512857e-12 * pft * pft ) * pfh & 742 - 0.164759e-06 * pfs & 743 +( 0.791325e-08 * pft - 0.933746e-06 ) * pft & 744 + 0.380374e-04 ) * pfh 745 746 ! 0. 2D fields of bottom temperature and salinity, and bottom slope 747 ! ----------------------------------------------------------------- 748 ! mbathy= number of w-level, minimum value=1 (cf domrea.F90) 749 # if defined key_vectopt_loop 750 jj = 1 751 DO ji = 1, jpij ! vector opt. (forced unrolling) 752 # else 753 DO jj = 1, jpj 754 DO ji = 1, jpi 755 # endif 756 ik = MAX( mbathy(ji,jj) - 1, 1 ) ! vertical index of the bottom ocean T-level 757 ztnb(ji,jj) = ptn(ji,jj,ik) * tmask(ji,jj,1) ! masked T and S at ocean bottom 758 zsnb(ji,jj) = psn(ji,jj,ik) * tmask(ji,jj,1) 759 zdep(ji,jj) = fsdept(ji,jj,ik) ! depth of the ocean bottom T-level 760 # if ! defined key_vectopt_loop 761 END DO 762 # endif 763 END DO 764 765 !!---------------------------------------------------------------------- 766 ! 1. Criteria of additional bottom diffusivity: grad(rho).grad(h)<0 767 ! -------------------------------------------- 768 ! Sign of the local density gradient along the i- and j-slopes 769 ! multiplied by the slope of the ocean bottom 770 771 SELECT CASE ( neos ) 772 773 CASE ( 0 ) ! Jackett and McDougall (1994) formulation 774 775 # if defined key_vectopt_loop 776 jj = 1 777 DO ji = 1, jpij-jpi ! vector opt. (forced unrolling) 778 # else 779 DO jj = 1, jpjm1 780 DO ji = 1, jpim1 781 # endif 782 ! temperature, salinity anomalie and depth 783 ztx = 0.5 * ( ztnb(ji,jj) + ztnb(ji+1,jj) ) 784 zsx = 0.5 * ( zsnb(ji,jj) + zsnb(ji+1,jj) ) - 35.0 785 zhx = 0.5 * ( zdep(ji,jj) + zdep(ji+1,jj) ) 786 ! 787 zty = 0.5 * ( ztnb(ji,jj+1) + ztnb(ji,jj) ) 788 zsy = 0.5 * ( zsnb(ji,jj+1) + zsnb(ji,jj) ) - 35.0 789 zhy = 0.5 * ( zdep(ji,jj+1) + zdep(ji,jj) ) 790 ! masked ratio alpha/beta 791 zalbetx = fsalbt( ztx, zsx, zhx ) * umask(ji,jj,1) 792 zalbety = fsalbt( zty, zsy, zhy ) * vmask(ji,jj,1) 793 ! local density gradient along i-bathymetric slope 794 zgdrhox = zalbetx * ( ztnb(ji+1,jj) - ztnb(ji,jj) ) & 795 - ( zsnb(ji+1,jj) - zsnb(ji,jj) ) 796 ! local density gradient along j-bathymetric slope 797 zgdrhoy = zalbety * ( ztnb(ji,jj+1) - ztnb(ji,jj) ) & 798 - ( zsnb(ji,jj+1) - zsnb(ji,jj) ) 799 ! sign of local i-gradient of density multiplied by the i-slope 800 pbblx(ji,jj) = 0.5 - SIGN( 0.5, -zgdrhox * ( zdep(ji+1,jj) - zdep(ji,jj) ) ) 801 ! sign of local j-gradient of density multiplied by the j-slope 802 pbbly(ji,jj) = 0.5 - SIGN( 0.5, -zgdrhoy * ( zdep(ji,jj+1) - zdep(ji,jj) ) ) 803 # if ! defined key_vectopt_loop 804 END DO 805 # endif 806 END DO 807 808 CASE ( 1 ) ! Linear formulation function of temperature only 809 ! 810 # if defined key_vectopt_loop 811 jj = 1 812 DO ji = 1, jpij-jpi ! vector opt. (forced unrolling) 813 # else 814 DO jj = 1, jpjm1 815 DO ji = 1, jpim1 816 # endif 817 ! local 'density/temperature' gradient along i-bathymetric slope 818 zgdrhox = ztnb(ji+1,jj) - ztnb(ji,jj) 819 ! local density gradient along j-bathymetric slope 820 zgdrhoy = ztnb(ji,jj+1) - ztnb(ji,jj) 821 ! sign of local i-gradient of density multiplied by the i-slope 822 pbblx(ji,jj) = 0.5 - SIGN( 0.5, -zgdrhox * ( zdep(ji+1,jj) - zdep(ji,jj) ) ) 823 ! sign of local j-gradient of density multiplied by the j-slope 824 pbbly(ji,jj) = 0.5 - SIGN( 0.5, -zgdrhoy * ( zdep(ji,jj+1) - zdep(ji,jj) ) ) 825 # if ! defined key_vectopt_loop 826 END DO 827 # endif 828 END DO 829 830 CASE ( 2 ) ! Linear formulation function of temperature and salinity 831 832 # if defined key_vectopt_loop 833 jj = 1 834 DO ji = 1, jpij-jpi ! vector opt. (forced unrolling) 835 # else 836 DO jj = 1, jpjm1 837 DO ji = 1, jpim1 838 # endif 839 ! local density gradient along i-bathymetric slope 840 zgdrhox = - ( rbeta*( zsnb(ji+1,jj) - zsnb(ji,jj) ) & 841 - ralpha*( ztnb(ji+1,jj) - ztnb(ji,jj) ) ) 842 ! local density gradient along j-bathymetric slope 843 zgdrhoy = - ( rbeta*( zsnb(ji,jj+1) - zsnb(ji,jj) ) & 844 - ralpha*( ztnb(ji,jj+1) - ztnb(ji,jj) ) ) 845 ! sign of local i-gradient of density multiplied by the i-slope 846 pbblx(ji,jj) = 0.5 - SIGN( 0.5, - zgdrhox * ( zdep(ji+1,jj) - zdep(ji,jj) ) ) 847 ! sign of local j-gradient of density multiplied by the j-slope 848 pbbly(ji,jj) = 0.5 - SIGN( 0.5, -zgdrhoy * ( zdep(ji,jj+1) - zdep(ji,jj) ) ) 849 # if ! defined key_vectopt_loop 850 END DO 851 # endif 852 END DO 853 854 CASE DEFAULT 855 856 WRITE(ctmp1,*) ' bad flag value for neos = ', neos 857 CALL ctl_stop(ctmp1) 858 859 END SELECT 860 861 ! Lateral boundary conditions 862 CALL lbc_lnk( pbblx, 'U', 1. ) 863 CALL lbc_lnk( pbbly, 'V', 1. ) 864 865 END SUBROUTINE bbl_sign 866 867 #endif 868 869 SUBROUTINE swap_dyn_data 870 !!---------------------------------------------------------------------- 871 !! *** ROUTINE swap_dyn_data *** 872 !! 873 !! ** Purpose : swap array data 874 !! 875 !! History : 876 !! 9.0 ! 07-09 (C. Ethe) 877 !!---------------------------------------------------------------------- 878 879 880 ! swap from record 2 to 1 881 tdta (:,:,:,1) = tdta (:,:,:,2) 882 sdta (:,:,:,1) = sdta (:,:,:,2) 883 avtdta (:,:,:,1) = avtdta (:,:,:,2) 884 udta (:,:,:,1) = udta (:,:,:,2) 885 vdta (:,:,:,1) = vdta (:,:,:,2) 886 wdta (:,:,:,1) = wdta (:,:,:,2) 887 #if defined key_trc_diatrd 888 hdivdta(:,:,:,1) = hdivdta(:,:,:,2) 889 #endif 890 891 #if defined key_ldfslp 892 uslpdta (:,:,:,1) = uslpdta (:,:,:,2) 893 vslpdta (:,:,:,1) = vslpdta (:,:,:,2) 894 wslpidta(:,:,:,1) = wslpidta(:,:,:,2) 895 wslpjdta(:,:,:,1) = wslpjdta(:,:,:,2) 896 #endif 897 hmlddta(:,:,1) = hmlddta(:,:,2) 898 wspddta(:,:,1) = wspddta(:,:,2) 899 frlddta(:,:,1) = frlddta(:,:,2) 900 empdta (:,:,1) = empdta (:,:,2) 901 qsrdta (:,:,1) = qsrdta (:,:,2) 902 903 #if ! defined key_off_degrad && defined key_traldf_c2d 904 ahtwdta(:,:,1) = ahtwdta(:,:,2) 905 # if defined key_trcldf_eiv 906 aeiwdta(:,:,1) = aeiwdta(:,:,2) 907 # endif 908 #endif 909 910 #if defined key_off_degrad 911 ahtudta(:,:,:,1) = ahtudta(:,:,:,2) 912 ahtvdta(:,:,:,1) = ahtvdta(:,:,:,2) 913 ahtwdta(:,:,:,1) = ahtwdta(:,:,:,2) 914 # if defined key_trcldf_eiv 915 aeiudta(:,:,:,1) = aeiudta(:,:,:,2) 916 aeivdta(:,:,:,1) = aeivdta(:,:,:,2) 917 aeiwdta(:,:,:,1) = aeiwdta(:,:,:,2) 918 # endif 919 #endif 920 921 #if defined key_trcbbl_dif || defined key_trcbbl_adv 922 bblxdta(:,:,1) = bblxdta(:,:,2) 923 bblydta(:,:,1) = bblydta(:,:,2) 924 #endif 925 926 END SUBROUTINE swap_dyn_data 927 928 SUBROUTINE assign_dyn_data 929 !!---------------------------------------------------------------------- 930 !! *** ROUTINE assign_dyn_data *** 931 !! 932 !! ** Purpose : Assign dynamical data to the data that have been read 933 !! without time interpolation 934 !! 935 !!---------------------------------------------------------------------- 936 937 tn (:,:,:) = tdta (:,:,:,2) 938 sn (:,:,:) = sdta (:,:,:,2) 939 avt(:,:,:) = avtdta(:,:,:,2) 940 941 un (:,:,:) = udta (:,:,:,2) 942 vn (:,:,:) = vdta (:,:,:,2) 943 wn (:,:,:) = wdta (:,:,:,2) 944 945 #if defined key_trc_diatrd 946 hdivn(:,:,:) = hdivdta(:,:,:,2) 947 #endif 948 949 #if defined key_zdfddm 950 avs(:,:,:) = avtdta (:,:,:,2) 951 #endif 952 953 954 #if defined key_ldfslp 955 uslp (:,:,:) = uslpdta (:,:,:,2) 956 vslp (:,:,:) = vslpdta (:,:,:,2) 957 wslpi(:,:,:) = wslpidta(:,:,:,2) 958 wslpj(:,:,:) = wslpjdta(:,:,:,2) 959 #endif 960 961 hmld(:,:) = hmlddta(:,:,2) 962 wndm(:,:) = wspddta(:,:,2) 963 fr_i(:,:) = frlddta(:,:,2) 964 emp (:,:) = empdta (:,:,2) 965 emps(:,:) = emp(:,:) 966 qsr (:,:) = qsrdta (:,:,2) 967 968 #if ! defined key_off_degrad && defined key_traldf_c2d 969 ahtw(:,:) = ahtwdta(:,:,2) 970 # if defined key_trcldf_eiv 971 aeiw(:,:) = aeiwdta(:,:,2) 972 # endif 973 #endif 974 975 #if defined key_off_degrad 976 ahtu(:,:,:) = ahtudta(:,:,:,2) 977 ahtv(:,:,:) = ahtvdta(:,:,:,2) 978 ahtw(:,:,:) = ahtwdta(:,:,:,2) 979 # if defined key_trcldf_eiv 980 aeiu(:,:,:) = aeiudta(:,:,:,2) 981 aeiv(:,:,:) = aeivdta(:,:,:,2) 982 aeiw(:,:,:) = aeiwdta(:,:,:,2) 983 # endif 984 985 #endif 986 987 #if defined key_trcbbl_dif || defined key_trcbbl_adv 988 bblx(:,:) = bblxdta(:,:,2) 989 bbly(:,:) = bblydta(:,:,2) 990 #endif 991 992 END SUBROUTINE assign_dyn_data 993 994 SUBROUTINE linear_interp_dyn_data( pweigh ) 995 !!---------------------------------------------------------------------- 996 !! *** ROUTINE linear_interp_dyn_data *** 997 !! 998 !! ** Purpose : linear interpolation of data 999 !! 1000 !!---------------------------------------------------------------------- 1001 !! * Argument 1002 REAL(wp), INTENT( in ) :: pweigh ! weigh 1003 1004 !! * Local declarations 1005 REAL(wp) :: zweighm1 1006 !!---------------------------------------------------------------------- 1007 1008 zweighm1 = 1. - pweigh 1009 1010 tn (:,:,:) = zweighm1 * tdta (:,:,:,1) + pweigh * tdta (:,:,:,2) 1011 sn (:,:,:) = zweighm1 * sdta (:,:,:,1) + pweigh * sdta (:,:,:,2) 1012 avt(:,:,:) = zweighm1 * avtdta(:,:,:,1) + pweigh * avtdta(:,:,:,2) 1013 1014 un (:,:,:) = zweighm1 * udta (:,:,:,1) + pweigh * udta (:,:,:,2) 1015 vn (:,:,:) = zweighm1 * vdta (:,:,:,1) + pweigh * vdta (:,:,:,2) 1016 wn (:,:,:) = zweighm1 * wdta (:,:,:,1) + pweigh * wdta (:,:,:,2) 1017 1018 #if defined key_trc_diatrd 1019 hdivn(:,:,:) = zweighm1 * hdivdta(:,:,:,1) + pweigh * hdivdta(:,:,:,2) 1020 #endif 1021 1022 #if defined key_zdfddm 1023 avs(:,:,:) = zweighm1 * avtdta (:,:,:,1) + pweigh * avtdta (:,:,:,2) 1024 #endif 1025 1026 1027 #if defined key_ldfslp 1028 uslp (:,:,:) = zweighm1 * uslpdta (:,:,:,1) + pweigh * uslpdta (:,:,:,2) 1029 vslp (:,:,:) = zweighm1 * vslpdta (:,:,:,1) + pweigh * vslpdta (:,:,:,2) 1030 wslpi(:,:,:) = zweighm1 * wslpidta(:,:,:,1) + pweigh * wslpidta(:,:,:,2) 1031 wslpj(:,:,:) = zweighm1 * wslpjdta(:,:,:,1) + pweigh * wslpjdta(:,:,:,2) 1032 #endif 1033 1034 hmld(:,:) = zweighm1 * hmlddta(:,:,1) + pweigh * hmlddta(:,:,2) 1035 wndm(:,:) = zweighm1 * wspddta(:,:,1) + pweigh * wspddta(:,:,2) 1036 fr_i(:,:) = zweighm1 * frlddta(:,:,1) + pweigh * frlddta(:,:,2) 1037 emp (:,:) = zweighm1 * empdta (:,:,1) + pweigh * empdta (:,:,2) 1038 emps(:,:) = emp(:,:) 1039 qsr (:,:) = zweighm1 * qsrdta (:,:,1) + pweigh * qsrdta (:,:,2) 1040 1041 #if ! defined key_off_degrad && defined key_traldf_c2d 1042 ahtw(:,:) = zweighm1 * ahtwdta(:,:,1) + pweigh * ahtwdta(:,:,2) 1043 # if defined key_trcldf_eiv 1044 aeiw(:,:) = zweighm1 * aeiwdta(:,:,1) + pweigh * aeiwdta(:,:,2) 1045 # endif 1046 #endif 1047 1048 #if defined key_off_degrad 1049 ahtu(:,:,:) = zweighm1 * ahtudta(:,:,:,1) + pweigh * ahtudta(:,:,:,2) 1050 ahtv(:,:,:) = zweighm1 * ahtvdta(:,:,:,1) + pweigh * ahtvdta(:,:,:,2) 1051 ahtw(:,:,:) = zweighm1 * ahtwdta(:,:,:,1) + pweigh * ahtwdta(:,:,:,2) 1052 # if defined key_trcldf_eiv 1053 aeiu(:,:,:) = zweighm1 * aeiudta(:,:,:,1) + pweigh * aeiudta(:,:,:,2) 1054 aeiv(:,:,:) = zweighm1 * aeivdta(:,:,:,1) + pweigh * aeivdta(:,:,:,2) 1055 aeiw(:,:,:) = zweighm1 * aeiwdta(:,:,:,1) + pweigh * aeiwdta(:,:,:,2) 1056 # endif 1057 #endif 1058 1059 #if defined key_trcbbl_dif || defined key_trcbbl_adv 1060 bblx(:,:) = zweighm1 * bblxdta(:,:,1) + pweigh * bblxdta(:,:,2) 1061 bbly(:,:) = zweighm1 * bblydta(:,:,1) + pweigh * bblydta(:,:,2) 1062 #endif 1063 1064 END SUBROUTINE linear_interp_dyn_data 908 1065 909 1066 END MODULE dtadyn
Note: See TracChangeset
for help on using the changeset viewer.