Changeset 7343 for branches/2015/dev_r5936_INGV1_WAVE/NEMOGCM/NEMO
- Timestamp:
- 2016-11-25T17:56:42+01:00 (8 years ago)
- Location:
- branches/2015/dev_r5936_INGV1_WAVE/NEMOGCM/NEMO/OPA_SRC
- Files:
-
- 5 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/2015/dev_r5936_INGV1_WAVE/NEMOGCM/NEMO/OPA_SRC/SBC/cpl_oasis3.F90
r5836 r7343 66 66 INTEGER :: nsnd ! total number of fields sent 67 67 INTEGER :: ncplmodel ! Maximum number of models to/from which NEMO is potentialy sending/receiving data 68 INTEGER, PUBLIC, PARAMETER :: nmaxfld=5 0! Maximum number of coupling fields68 INTEGER, PUBLIC, PARAMETER :: nmaxfld=55 ! Maximum number of coupling fields 69 69 INTEGER, PUBLIC, PARAMETER :: nmaxcat=5 ! Maximum number of coupling fields 70 70 INTEGER, PUBLIC, PARAMETER :: nmaxcpl=5 ! Maximum number of coupling fields -
branches/2015/dev_r5936_INGV1_WAVE/NEMOGCM/NEMO/OPA_SRC/SBC/sbccpl.F90
r5836 r7343 23 23 USE sbcapr 24 24 USE sbcdcy ! surface boundary condition: diurnal cycle 25 USE sbcwave ! surface boundary condition: waves 25 26 USE phycst ! physical constants 26 27 #if defined key_lim3 … … 105 106 INTEGER, PARAMETER :: jpr_e3t1st = 41 ! first T level thickness 106 107 INTEGER, PARAMETER :: jpr_fraqsr = 42 ! fraction of solar net radiation absorbed in the first ocean level 107 INTEGER, PARAMETER :: jprcv = 42 ! total number of fields received 108 INTEGER, PARAMETER :: jpr_mslp = 43 ! mean sea level pressure 109 INTEGER, PARAMETER :: jpr_hsig = 44 ! Hsig 110 INTEGER, PARAMETER :: jpr_phioc = 45 ! Wave=>ocean energy flux 111 INTEGER, PARAMETER :: jpr_sdrftx = 46 ! Stokes drift on grid 1 112 INTEGER, PARAMETER :: jpr_sdrfty = 47 ! Stokes drift on grid 2 113 INTEGER, PARAMETER :: jpr_wper = 48 ! Mean wave period 114 INTEGER, PARAMETER :: jpr_wnum = 49 ! Mean wavenumber 115 INTEGER, PARAMETER :: jpr_wstrf = 50 ! Stress fraction adsorbed by waves 116 INTEGER, PARAMETER :: jpr_wdrag = 51 ! Neutral surface drag coefficient 117 INTEGER, PARAMETER :: jprcv = 51 ! total number of fields received 108 118 109 119 INTEGER, PARAMETER :: jps_fice = 1 ! ice fraction sent to the atmosphere … … 135 145 INTEGER, PARAMETER :: jps_e3t1st = 27 ! first level depth (vvl) 136 146 INTEGER, PARAMETER :: jps_fraqsr = 28 ! fraction of solar net radiation absorbed in the first ocean level 137 INTEGER, PARAMETER :: jpsnd = 28 ! total number of fields sended 147 INTEGER, PARAMETER :: jps_ficet = 29 ! total ice fraction 148 INTEGER, PARAMETER :: jps_ocxw = 30 ! currents on grid 1 149 INTEGER, PARAMETER :: jps_ocyw = 31 ! currents on grid 2 150 INTEGER, PARAMETER :: jps_wlev = 32 ! water level 151 INTEGER, PARAMETER :: jpsnd = 32 ! total number of fields sent 138 152 139 153 ! !!** namelist namsbc_cpl ** … … 149 163 ! Received from the atmosphere ! 150 164 TYPE(FLD_C) :: sn_rcv_w10m, sn_rcv_taumod, sn_rcv_tau, sn_rcv_dqnsdt, sn_rcv_qsr, sn_rcv_qns, sn_rcv_emp, sn_rcv_rnf 151 TYPE(FLD_C) :: sn_rcv_cal, sn_rcv_iceflx, sn_rcv_co2 165 TYPE(FLD_C) :: sn_rcv_cal, sn_rcv_iceflx, sn_rcv_co2, sn_rcv_mslp 166 ! Send to waves 167 TYPE(FLD_C) :: sn_snd_ifrac, sn_snd_crtw, sn_snd_wlev 168 ! Received from waves 169 TYPE(FLD_C) :: sn_rcv_hsig,sn_rcv_phioc,sn_rcv_sdrfx,sn_rcv_sdrfy,sn_rcv_wper,sn_rcv_wnum,sn_rcv_wstrf,sn_rcv_wdrag 152 170 ! Other namelist parameters ! 153 171 INTEGER :: nn_cplmodel ! Maximum number of models to/from which NEMO is potentialy sending/receiving data … … 161 179 162 180 REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: albedo_oce_mix ! ocean albedo sent to atmosphere (mix clear/overcast sky) 181 182 REAL(wp) :: rpref = 101000._wp ! reference atmospheric pressure[N/m2] 183 REAL(wp) :: r1_grau ! = 1.e0 / (grav * rau0) 163 184 164 185 INTEGER , ALLOCATABLE, SAVE, DIMENSION( :) :: nrcvinfo ! OASIS info argument … … 179 200 !! *** FUNCTION sbc_cpl_alloc *** 180 201 !!---------------------------------------------------------------------- 181 INTEGER :: ierr( 3)202 INTEGER :: ierr(4) 182 203 !!---------------------------------------------------------------------- 183 204 ierr(:) = 0 … … 190 211 ALLOCATE( xcplmask(jpi,jpj,0:nn_cplmodel) , STAT=ierr(3) ) 191 212 ! 213 IF( .NOT. ln_apr_dyn ) ALLOCATE( ssh_ib(jpi,jpj), ssh_ibb(jpi,jpj), apr(jpi, jpj), STAT=ierr(4) ) 214 192 215 sbc_cpl_alloc = MAXVAL( ierr ) 193 216 IF( lk_mpp ) CALL mpp_sum ( sbc_cpl_alloc ) … … 216 239 REAL(wp), POINTER, DIMENSION(:,:) :: zacs, zaos 217 240 !! 218 NAMELIST/namsbc_cpl/ sn_snd_temp, sn_snd_alb , sn_snd_thick, sn_snd_crt , sn_snd_co2, & 219 & sn_rcv_w10m, sn_rcv_taumod, sn_rcv_tau , sn_rcv_dqnsdt, sn_rcv_qsr, & 220 & sn_rcv_qns , sn_rcv_emp , sn_rcv_rnf , sn_rcv_cal , sn_rcv_iceflx, & 221 & sn_rcv_co2 , nn_cplmodel , ln_usecplmask 241 NAMELIST/namsbc_cpl/ sn_snd_temp , sn_snd_alb , sn_snd_thick , sn_snd_crt , sn_snd_co2, & 242 & sn_rcv_w10m, sn_rcv_taumod, sn_rcv_tau , sn_rcv_dqnsdt, sn_rcv_qsr, & 243 & sn_snd_ifrac, sn_snd_crtw , sn_snd_wlev , sn_rcv_hsig , sn_rcv_phioc , & 244 & sn_rcv_sdrfx, sn_rcv_sdrfy, sn_rcv_wper , sn_rcv_wnum , sn_rcv_wstrf , & 245 & sn_rcv_wdrag, sn_rcv_qns , sn_rcv_emp , sn_rcv_rnf , sn_rcv_cal , & 246 & sn_rcv_iceflx,sn_rcv_co2 , nn_cplmodel , ln_usecplmask, sn_rcv_mslp 222 247 !!--------------------------------------------------------------------- 223 248 ! … … 260 285 WRITE(numout,*)' sea ice heat fluxes = ', TRIM(sn_rcv_iceflx%cldes), ' (', TRIM(sn_rcv_iceflx%clcat), ')' 261 286 WRITE(numout,*)' atm co2 = ', TRIM(sn_rcv_co2%cldes ), ' (', TRIM(sn_rcv_co2%clcat ), ')' 287 WRITE(numout,*)' significant wave heigth = ', TRIM(sn_rcv_hsig%cldes ), ' (', TRIM(sn_rcv_hsig%clcat ), ')' 288 WRITE(numout,*)' wave to oce energy flux = ', TRIM(sn_rcv_phioc%cldes ), ' (', TRIM(sn_rcv_phioc%clcat ), ')' 289 WRITE(numout,*)' Surface Stokes drift grid u = ', TRIM(sn_rcv_sdrfx%cldes ), ' (', TRIM(sn_rcv_sdrfx%clcat ), ')' 290 WRITE(numout,*)' Surface Stokes drift grid v = ', TRIM(sn_rcv_sdrfy%cldes ), ' (', TRIM(sn_rcv_sdrfy%clcat ), ')' 291 WRITE(numout,*)' Mean wave period = ', TRIM(sn_rcv_wper%cldes ), ' (', TRIM(sn_rcv_wper%clcat ), ')' 292 WRITE(numout,*)' Mean wave number = ', TRIM(sn_rcv_wnum%cldes ), ' (', TRIM(sn_rcv_wnum%clcat ), ')' 293 WRITE(numout,*)' Stress frac adsorbed by waves = ', TRIM(sn_rcv_wstrf%cldes ), ' (', TRIM(sn_rcv_wstrf%clcat ), ')' 294 WRITE(numout,*)' Neutral surf drag coefficient = ', TRIM(sn_rcv_wdrag%cldes ), ' (', TRIM(sn_rcv_wdrag%clcat ), ')' 262 295 WRITE(numout,*)' sent fields (multiple ice categories)' 263 296 WRITE(numout,*)' surface temperature = ', TRIM(sn_snd_temp%cldes ), ' (', TRIM(sn_snd_temp%clcat ), ')' 264 297 WRITE(numout,*)' albedo = ', TRIM(sn_snd_alb%cldes ), ' (', TRIM(sn_snd_alb%clcat ), ')' 265 298 WRITE(numout,*)' ice/snow thickness = ', TRIM(sn_snd_thick%cldes ), ' (', TRIM(sn_snd_thick%clcat ), ')' 299 WRITE(numout,*)' total ice fraction = ', TRIM(sn_snd_ifrac%cldes ), ' (', TRIM(sn_snd_ifrac%clcat ), ')' 266 300 WRITE(numout,*)' surface current = ', TRIM(sn_snd_crt%cldes ), ' (', TRIM(sn_snd_crt%clcat ), ')' 267 301 WRITE(numout,*)' - referential = ', sn_snd_crt%clvref … … 269 303 WRITE(numout,*)' - mesh = ', sn_snd_crt%clvgrd 270 304 WRITE(numout,*)' oce co2 flux = ', TRIM(sn_snd_co2%cldes ), ' (', TRIM(sn_snd_co2%clcat ), ')' 305 WRITE(numout,*)' water level = ', TRIM(sn_snd_wlev%cldes ), ' (', TRIM(sn_snd_wlev%clcat ), ')' 306 WRITE(numout,*)' mean sea level pressure = ', TRIM(sn_rcv_mslp%cldes ), ' (', TRIM(sn_rcv_mslp%clcat ), ')' 307 WRITE(numout,*)' surface current to waves = ', TRIM(sn_snd_crtw%cldes ), ' (', TRIM(sn_snd_crtw%clcat ), ')' 308 WRITE(numout,*)' - referential = ', sn_snd_crtw%clvref 309 WRITE(numout,*)' - orientation = ', sn_snd_crtw%clvor 310 WRITE(numout,*)' - mesh = ', sn_snd_crtw%clvgrd 271 311 WRITE(numout,*)' nn_cplmodel = ', nn_cplmodel 272 312 WRITE(numout,*)' ln_usecplmask = ', ln_usecplmask … … 307 347 ! 308 348 ! Vectors: change of sign at north fold ONLY if on the local grid 349 IF( TRIM( sn_rcv_tau%cldes ) == 'oce only' .OR. TRIM(sn_rcv_tau%cldes ) == 'oce and ice') THEN ! avoid working with the atmospheric fields if they are not coupled 309 350 IF( TRIM( sn_rcv_tau%clvor ) == 'local grid' ) srcv(jpr_otx1:jpr_itz2)%nsgn = -1. 310 351 … … 374 415 srcv(jpr_ity1)%clgrid = 'V' ! i.e. it is always at U- & V-points for i- & j-comp. resp. 375 416 ENDIF 417 ENDIF 376 418 377 419 ! ! ------------------------- ! … … 470 512 ! ! ------------------------- ! 471 513 srcv(jpr_co2 )%clname = 'O_AtmCO2' ; IF( TRIM(sn_rcv_co2%cldes ) == 'coupled' ) srcv(jpr_co2 )%laction = .TRUE. 514 515 ! ! ------------------------- ! 516 ! ! Mean Sea Level Pressure ! 517 ! ! ------------------------- ! 518 srcv(jpr_mslp)%clname = 'O_MSLP' ; IF( TRIM(sn_rcv_mslp%cldes ) == 'coupled' ) srcv(jpr_mslp)%laction = .TRUE. 519 472 520 ! ! ------------------------- ! 473 521 ! ! topmelt and botmelt ! … … 483 531 srcv(jpr_topm:jpr_botm)%laction = .TRUE. 484 532 ENDIF 533 ! ! ------------------------- ! 534 ! ! Wave breaking ! 535 ! ! ------------------------- ! 536 srcv(jpr_hsig)%clname = 'O_Hsigwa' ! significant wave height 537 IF( TRIM(sn_rcv_hsig%cldes ) == 'coupled' ) THEN 538 srcv(jpr_hsig)%laction = .TRUE. 539 cpl_hsig = .TRUE. 540 ENDIF 541 srcv(jpr_phioc)%clname = 'O_PhiOce' ! wave to ocean energy 542 IF( TRIM(sn_rcv_phioc%cldes ) == 'coupled' ) THEN 543 srcv(jpr_phioc)%laction = .TRUE. 544 cpl_phioc = .TRUE. 545 ENDIF 546 srcv(jpr_sdrftx)%clname = 'O_Sdrfx' ! Stokes drift in the u direction 547 IF( TRIM(sn_rcv_sdrfx%cldes ) == 'coupled' ) THEN 548 srcv(jpr_sdrftx)%laction = .TRUE. 549 cpl_sdrftx = .TRUE. 550 ENDIF 551 srcv(jpr_sdrfty)%clname = 'O_Sdrfy' ! Stokes drift in the v direction 552 IF( TRIM(sn_rcv_sdrfy%cldes ) == 'coupled' ) THEN 553 srcv(jpr_sdrfty)%laction = .TRUE. 554 cpl_sdrfty = .TRUE. 555 ENDIF 556 srcv(jpr_wper)%clname = 'O_WPer' ! mean wave period 557 IF( TRIM(sn_rcv_wper%cldes ) == 'coupled' ) THEN 558 srcv(jpr_wper)%laction = .TRUE. 559 cpl_wper = .TRUE. 560 ENDIF 561 srcv(jpr_wnum)%clname = 'O_WNum' ! mean wave number 562 IF( TRIM(sn_rcv_wnum%cldes ) == 'coupled' ) THEN 563 srcv(jpr_wnum)%laction = .TRUE. 564 cpl_wnum = .TRUE. 565 ENDIF 566 srcv(jpr_wstrf)%clname = 'O_WStrf' ! stress fraction adsorbed by the wave 567 IF( TRIM(sn_rcv_wstrf%cldes ) == 'coupled' ) THEN 568 srcv(jpr_wstrf)%laction = .TRUE. 569 cpl_wstrf = .TRUE. 570 ENDIF 571 srcv(jpr_wdrag)%clname = 'O_WDrag' ! neutral surface drag coefficient 572 IF( TRIM(sn_rcv_wdrag%cldes ) == 'coupled' ) THEN 573 srcv(jpr_wdrag)%laction = .TRUE. 574 cpl_wdrag = .TRUE. 575 ENDIF 576 ! 485 577 ! ! ------------------------------- ! 486 578 ! ! OPA-SAS coupling - rcv by opa ! … … 637 729 ! ! ------------------------- ! 638 730 ssnd(jps_fice)%clname = 'OIceFrc' 731 ssnd(jps_ficet)%clname = 'OIceFrcT' 639 732 ssnd(jps_hice)%clname = 'OIceTck' 640 733 ssnd(jps_hsnw)%clname = 'OSnwTck' … … 645 738 ENDIF 646 739 740 IF (TRIM( sn_snd_ifrac%cldes ) == 'coupled') ssnd(jps_ficet)%laction = .TRUE. 741 647 742 SELECT CASE ( TRIM( sn_snd_thick%cldes ) ) 648 743 CASE( 'none' ) ! nothing to do … … 665 760 ssnd(jps_ocy1)%clname = 'O_OCury1' ; ssnd(jps_ivy1)%clname = 'O_IVely1' 666 761 ssnd(jps_ocz1)%clname = 'O_OCurz1' ; ssnd(jps_ivz1)%clname = 'O_IVelz1' 762 ssnd(jps_ocxw)%clname = 'O_OCurxw' 763 ssnd(jps_ocyw)%clname = 'O_OCuryw' 667 764 ! 668 765 ssnd(jps_ocx1:jps_ivz1)%nsgn = -1. ! vectors: change of the sign at the north fold … … 685 782 END SELECT 686 783 784 ssnd(jps_ocxw:jps_ocyw)%nsgn = -1. ! vectors: change of the sign at the north fold 785 786 IF( sn_snd_crtw%clvgrd == 'U,V' ) THEN 787 ssnd(jps_ocxw)%clgrid = 'U' ; ssnd(jps_ocyw)%clgrid = 'V' 788 ELSE IF( sn_snd_crtw%clvgrd /= 'T' ) THEN 789 CALL ctl_stop( 'sn_snd_crtw%clvgrd must be equal to T' ) 790 ENDIF 791 IF( TRIM( sn_snd_crtw%clvor ) == 'eastward-northward' ) ssnd(jps_ocxw:jps_ocyw)%nsgn = 1. 792 SELECT CASE( TRIM( sn_snd_crtw%cldes ) ) 793 CASE( 'none' ) ; ssnd(jps_ocxw:jps_ocyw)%laction = .FALSE. 794 CASE( 'oce only' ) ; ssnd(jps_ocxw:jps_ocyw)%laction = .TRUE. 795 CASE( 'weighted oce and ice' ) ! nothing to do 796 CASE( 'mixed oce-ice' ) ; ssnd(jps_ivx1:jps_ivz1)%laction = .FALSE. 797 CASE default ; CALL ctl_stop( 'sbc_cpl_init: wrong definition of sn_snd_crtw%cldes' ) 798 END SELECT 799 687 800 ! ! ------------------------- ! 688 801 ! ! CO2 flux ! 689 802 ! ! ------------------------- ! 690 803 ssnd(jps_co2)%clname = 'O_CO2FLX' ; IF( TRIM(sn_snd_co2%cldes) == 'coupled' ) ssnd(jps_co2 )%laction = .TRUE. 804 805 ! ! ------------------------- ! 806 ! ! Sea surface height ! 807 ! ! ------------------------- ! 808 ssnd(jps_wlev)%clname = 'O_Wlevel' ; IF( TRIM(sn_snd_wlev%cldes) == 'coupled' ) ssnd(jps_wlev)%laction = .TRUE. 691 809 692 810 ! ! ------------------------------- ! … … 783 901 IF( ln_dm2dc .AND. ln_cpl .AND. ncpl_qsr_freq /= 86400 ) & 784 902 & CALL ctl_stop( 'sbc_cpl_init: diurnal cycle reconstruction (ln_dm2dc) needs daily couping for solar radiation' ) 785 ncpl_qsr_freq = 86400 / ncpl_qsr_freq903 IF( ln_dm2dc .AND. ln_cpl ) ncpl_qsr_freq = 86400 / ncpl_qsr_freq 786 904 787 905 CALL wrk_dealloc( jpi,jpj, zacs, zaos ) … … 837 955 !! emp upward mass flux [evap. - precip. (- runoffs) (- calving)] (ocean only case) 838 956 !!---------------------------------------------------------------------- 957 USE zdf_oce, ONLY : ln_zdfqiao 958 959 IMPLICIT NONE 960 839 961 INTEGER, INTENT(in) :: kt ! ocean model time step index 840 962 INTEGER, INTENT(in) :: k_fsbc ! frequency of sbc (-> ice model) computation … … 992 1114 IF( srcv(jpr_co2)%laction ) atm_co2(:,:) = frcv(jpr_co2)%z3(:,:,1) 993 1115 #endif 1116 ! 1117 ! ! ========================= ! 1118 ! ! Mean Sea Level Pressure ! (taum) 1119 ! ! ========================= ! 1120 ! 1121 IF( srcv(jpr_mslp)%laction ) THEN ! UKMO SHELF effect of atmospheric pressure on SSH 1122 IF( kt /= nit000 ) ssh_ibb(:,:) = ssh_ib(:,:) !* Swap of ssh_ib fields 1123 1124 r1_grau = 1.e0 / (grav * rau0) !* constant for optimization 1125 ssh_ib(:,:) = - ( frcv(jpr_mslp)%z3(:,:,1) - rpref ) * r1_grau ! equivalent ssh (inverse barometer) 1126 apr (:,:) = frcv(jpr_mslp)%z3(:,:,1) !atmospheric pressure 1127 1128 IF( kt == nit000 ) ssh_ibb(:,:) = ssh_ib(:,:) ! correct this later (read from restart if possible) 1129 END IF 1130 ! 1131 IF( ln_sdw ) THEN ! Stokes Drift correction activated 1132 ! ! ========================= ! 1133 ! ! Stokes drift u ! 1134 ! ! ========================= ! 1135 IF( srcv(jpr_sdrftx)%laction ) zusd2dt(:,:) = frcv(jpr_sdrftx)%z3(:,:,1) 1136 ! 1137 ! ! ========================= ! 1138 ! ! Stokes drift v ! 1139 ! ! ========================= ! 1140 IF( srcv(jpr_sdrfty)%laction ) zvsd2dt(:,:) = frcv(jpr_sdrfty)%z3(:,:,1) 1141 ! 1142 ! ! ========================= ! 1143 ! ! Wave mean period ! 1144 ! ! ========================= ! 1145 IF( srcv(jpr_wper)%laction ) wmp(:,:) = frcv(jpr_wper)%z3(:,:,1) 1146 ! 1147 ! ! ========================= ! 1148 ! ! Significant wave height ! 1149 ! ! ========================= ! 1150 IF( srcv(jpr_hsig)%laction ) swh(:,:) = frcv(jpr_hsig)%z3(:,:,1) 1151 ! 1152 ! ! ========================= ! 1153 ! ! Vertical mixing Qiao ! 1154 ! ! ========================= ! 1155 IF( srcv(jpr_wnum)%laction .AND. ln_zdfqiao ) wnum(:,:) = frcv(jpr_wnum)%z3(:,:,1) 1156 1157 ! Calculate the 3D Stokes drift both in coupled and not fully uncoupled mode 1158 IF( srcv(jpr_sdrftx)%laction .OR. srcv(jpr_sdrfty)%laction .OR. srcv(jpr_wper)%laction & 1159 .OR. srcv(jpr_hsig)%laction ) THEN 1160 CALL sbc_stokes() 1161 IF( ln_zdfqiao .AND. .NOT. srcv(jpr_wnum)%laction ) CALL sbc_qiao() 1162 ENDIF 1163 IF( ln_zdfqiao .AND. srcv(jpr_wnum)%laction ) CALL sbc_qiao() 1164 ENDIF 1165 ! ! ========================= ! 1166 ! ! Stress adsorbed by waves ! 1167 ! ! ========================= ! 1168 IF( srcv(jpr_wstrf)%laction .AND. ln_tauoc ) tauoc_wave(:,:) = frcv(jpr_wstrf)%z3(:,:,1) 1169 1170 ! ! ========================= ! 1171 ! ! Wave drag coefficient ! 1172 ! ! ========================= ! 1173 IF( srcv(jpr_wdrag)%laction .AND. ln_cdgw ) cdn_wave(:,:) = frcv(jpr_wdrag)%z3(:,:,1) 994 1174 995 1175 ! Fields received by SAS when OASIS coupling … … 1984 2164 ENDIF 1985 2165 ! 2166 ! ! ------------------------- ! 2167 ! ! Surface current to waves ! 2168 ! ! ------------------------- ! 2169 IF( ssnd(jps_ocxw)%laction .OR. ssnd(jps_ocyw)%laction ) THEN 2170 ! 2171 ! j+1 j -----V---F 2172 ! surface velocity always sent from T point ! | 2173 ! j | T U 2174 ! | | 2175 ! j j-1 -I-------| 2176 ! (for I) | | 2177 ! i-1 i i 2178 ! i i+1 (for I) 2179 SELECT CASE( TRIM( sn_snd_crtw%cldes ) ) 2180 CASE( 'oce only' ) ! C-grid ==> T 2181 DO jj = 2, jpjm1 2182 DO ji = fs_2, fs_jpim1 ! vector opt. 2183 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) 2184 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji , jj-1,1) ) 2185 END DO 2186 END DO 2187 CASE( 'weighted oce and ice' ) 2188 SELECT CASE ( cp_ice_msh ) 2189 CASE( 'C' ) ! Ocean and Ice on C-grid ==> T 2190 DO jj = 2, jpjm1 2191 DO ji = fs_2, fs_jpim1 ! vector opt. 2192 zotx1(ji,jj) = 0.5 * ( un (ji,jj,1) + un (ji-1,jj ,1) ) * zfr_l(ji,jj) 2193 zoty1(ji,jj) = 0.5 * ( vn (ji,jj,1) + vn (ji ,jj-1,1) ) * zfr_l(ji,jj) 2194 zitx1(ji,jj) = 0.5 * ( u_ice(ji,jj ) + u_ice(ji-1,jj ) ) * fr_i(ji,jj) 2195 zity1(ji,jj) = 0.5 * ( v_ice(ji,jj ) + v_ice(ji ,jj-1 ) ) * fr_i(ji,jj) 2196 END DO 2197 END DO 2198 CASE( 'I' ) ! Ocean on C grid, Ice on I-point (B-grid) ==> T 2199 DO jj = 2, jpjm1 2200 DO ji = 2, jpim1 ! NO vector opt. 2201 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) 2202 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) 2203 zitx1(ji,jj) = 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1) & 2204 & + u_ice(ji+1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 2205 zity1(ji,jj) = 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1) & 2206 & + v_ice(ji+1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 2207 END DO 2208 END DO 2209 CASE( 'F' ) ! Ocean on C grid, Ice on F-point (B-grid) ==> T 2210 DO jj = 2, jpjm1 2211 DO ji = 2, jpim1 ! NO vector opt. 2212 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) 2213 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) 2214 zitx1(ji,jj) = 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1) & 2215 & + u_ice(ji-1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 2216 zity1(ji,jj) = 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1) & 2217 & + v_ice(ji-1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 2218 END DO 2219 END DO 2220 END SELECT 2221 CALL lbc_lnk( zitx1, 'T', -1. ) ; CALL lbc_lnk( zity1, 'T', -1. ) 2222 CASE( 'mixed oce-ice' ) 2223 SELECT CASE ( cp_ice_msh ) 2224 CASE( 'C' ) ! Ocean and Ice on C-grid ==> T 2225 DO jj = 2, jpjm1 2226 DO ji = fs_2, fs_jpim1 ! vector opt. 2227 zotx1(ji,jj) = 0.5 * ( un (ji,jj,1) + un (ji-1,jj ,1) ) * zfr_l(ji,jj) & 2228 & + 0.5 * ( u_ice(ji,jj ) + u_ice(ji-1,jj ) ) * fr_i(ji,jj) 2229 zoty1(ji,jj) = 0.5 * ( vn (ji,jj,1) + vn (ji ,jj-1,1) ) * zfr_l(ji,jj) & 2230 & + 0.5 * ( v_ice(ji,jj ) + v_ice(ji ,jj-1 ) ) * fr_i(ji,jj) 2231 END DO 2232 END DO 2233 CASE( 'I' ) ! Ocean on C grid, Ice on I-point (B-grid) ==> T 2234 DO jj = 2, jpjm1 2235 DO ji = 2, jpim1 ! NO vector opt. 2236 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) & 2237 & + 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1) & 2238 & + u_ice(ji+1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 2239 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) & 2240 & + 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1) & 2241 & + v_ice(ji+1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 2242 END DO 2243 END DO 2244 CASE( 'F' ) ! Ocean on C grid, Ice on F-point (B-grid) ==> T 2245 DO jj = 2, jpjm1 2246 DO ji = 2, jpim1 ! NO vector opt. 2247 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) & 2248 & + 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1) & 2249 & + u_ice(ji-1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 2250 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) & 2251 & + 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1) & 2252 & + v_ice(ji-1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 2253 END DO 2254 END DO 2255 END SELECT 2256 END SELECT 2257 CALL lbc_lnk( zotx1, ssnd(jps_ocxw)%clgrid, -1. ) ; CALL lbc_lnk( zoty1, ssnd(jps_ocyw)%clgrid, -1. ) 2258 ! 2259 ! 2260 IF( TRIM( sn_snd_crtw%clvor ) == 'eastward-northward' ) THEN ! Rotation of the components 2261 ! ! Ocean component 2262 CALL rot_rep( zotx1, zoty1, ssnd(jps_ocxw)%clgrid, 'ij->e', ztmp1 ) ! 1st component 2263 CALL rot_rep( zotx1, zoty1, ssnd(jps_ocxw)%clgrid, 'ij->n', ztmp2 ) ! 2nd component 2264 zotx1(:,:) = ztmp1(:,:) ! overwrite the components 2265 zoty1(:,:) = ztmp2(:,:) 2266 IF( ssnd(jps_ivx1)%laction ) THEN ! Ice component 2267 CALL rot_rep( zitx1, zity1, ssnd(jps_ivx1)%clgrid, 'ij->e', ztmp1 ) ! 1st component 2268 CALL rot_rep( zitx1, zity1, ssnd(jps_ivx1)%clgrid, 'ij->n', ztmp2 ) ! 2nd component 2269 zitx1(:,:) = ztmp1(:,:) ! overwrite the components 2270 zity1(:,:) = ztmp2(:,:) 2271 ENDIF 2272 ENDIF 2273 ! 2274 ! ! spherical coordinates to cartesian -> 2 components to 3 components 2275 ! IF( TRIM( sn_snd_crtw%clvref ) == 'cartesian' ) THEN 2276 ! ztmp1(:,:) = zotx1(:,:) ! ocean currents 2277 ! ztmp2(:,:) = zoty1(:,:) 2278 ! CALL oce2geo ( ztmp1, ztmp2, 'T', zotx1, zoty1, zotz1 ) 2279 ! ! 2280 ! IF( ssnd(jps_ivx1)%laction ) THEN ! ice velocities 2281 ! ztmp1(:,:) = zitx1(:,:) 2282 ! ztmp1(:,:) = zity1(:,:) 2283 ! CALL oce2geo ( ztmp1, ztmp2, 'T', zitx1, zity1, zitz1 ) 2284 ! ENDIF 2285 ! ENDIF 2286 ! 2287 IF( ssnd(jps_ocxw)%laction ) CALL cpl_snd( jps_ocxw, isec, RESHAPE ( zotx1, (/jpi,jpj,1/) ), info ) ! ocean x current 1st grid 2288 IF( ssnd(jps_ocyw)%laction ) CALL cpl_snd( jps_ocyw, isec, RESHAPE ( zoty1, (/jpi,jpj,1/) ), info ) ! ocean y current 1st grid 2289 ! 2290 ENDIF 2291 ! 2292 IF( ssnd(jps_ficet)%laction ) THEN 2293 CALL cpl_snd( jps_ficet, isec, RESHAPE ( fr_i, (/jpi,jpj,1/) ), info ) 2294 END IF 2295 ! ! ------------------------- ! 2296 ! ! Water levels to waves ! 2297 ! ! ------------------------- ! 2298 IF( ssnd(jps_wlev)%laction ) THEN 2299 IF( ln_apr_dyn ) THEN 2300 IF( kt /= nit000 ) THEN 2301 ztmp1(:,:) = sshb(:,:) - 0.5 * ( ssh_ib(:,:) + ssh_ibb(:,:) ) 2302 ELSE 2303 ztmp1(:,:) = sshb(:,:) 2304 ENDIF 2305 ELSE 2306 ztmp1(:,:) = sshn(:,:) 2307 ENDIF 2308 CALL cpl_snd( jps_wlev , isec, RESHAPE ( ztmp1, (/jpi,jpj,1/) ), info ) 2309 END IF 1986 2310 ! 1987 2311 ! Fields sent by OPA to SAS when doing OPA<->SAS coupling -
branches/2015/dev_r5936_INGV1_WAVE/NEMOGCM/NEMO/OPA_SRC/SBC/sbcmod.F90
r7221 r7343 131 131 WRITE(numout,*) ' ocean-atmosphere coupled formulation ln_cpl = ', ln_cpl 132 132 WRITE(numout,*) ' forced-coupled mixed formulation ln_mixcpl = ', ln_mixcpl 133 WRITE(numout,*) ' wave physics ln_wave = ', ln_wave 134 WRITE(numout,*) ' Stokes drift corr. to vert. velocity ln_sdw = ', ln_sdw 135 WRITE(numout,*) ' wave modified ocean stress ln_tauoc = ', ln_tauoc 136 WRITE(numout,*) ' Stokes coriolis term ln_stcor = ', ln_stcor 137 WRITE(numout,*) ' neutral drag coefficient (CORE, MFS) ln_cdgw = ', ln_cdgw 133 138 WRITE(numout,*) ' OASIS coupling (with atm or sas) lk_oasis = ', lk_oasis 134 139 WRITE(numout,*) ' components of your executable nn_components = ', nn_components … … 372 377 SELECT CASE( nsbc ) 373 378 CASE( 0,1,2,3,5,-1 ) ; 374 IF(lwp ) WRITE(numout,*) 'WARNING: You are subtracting the wave stress to the ocean. &379 IF(lwp .AND. kt == nit000 ) WRITE(numout,*) 'WARNING: You are subtracting the wave stress to the ocean. & 375 380 & If not requested select ln_tauoc=.false' 376 381 END SELECT -
branches/2015/dev_r5936_INGV1_WAVE/NEMOGCM/NEMO/OPA_SRC/SBC/sbcwave.F90
r7340 r7343 7 7 !! : 3.4 ! 2012-10 (M. Adani) Stokes Drift 8 8 !! 3.6 ! 2014-09 (E. Clementi,P. Oddo) New Stokes Drift Computation 9 !! 3.6 ! 2014-09 (Clementi E, Oddo P)New Stokes Drift Computation 9 10 !!---------------------------------------------------------------------- 10 11 … … 13 14 !!---------------------------------------------------------------------- 14 15 USE oce ! 15 USE sbc_oce 16 USE sbc_oce ! Surface boundary condition: ocean fields 16 17 USE bdy_oce ! 17 18 USE domvvl ! … … 19 20 USE in_out_manager ! I/O manager 20 21 USE lib_mpp ! distribued memory computing library 21 USE fldread 22 USE fldread ! read input fields 22 23 USE wrk_nemo ! 23 24 USE phycst ! physical constants … … 26 27 PRIVATE 27 28 29 PUBLIC sbc_stokes, sbc_qiao ! routines called in sbccpl 28 30 PUBLIC sbc_wave ! routine called in sbcmod 29 31 30 INTEGER, PARAMETER :: jpfld = 4 ! number of files to read for stokes drift 31 INTEGER, PARAMETER :: jp_usd = 1 ! index of stokes drift (i-component) (m/s) at T-point 32 INTEGER, PARAMETER :: jp_vsd = 2 ! index of stokes drift (j-component) (m/s) at T-point 33 INTEGER, PARAMETER :: jp_swh = 3 ! index of significant wave hight (m) at T-point 34 INTEGER, PARAMETER :: jp_wmp = 4 ! index of mean wave period (s) at T-point 32 ! Variables checking if the wave parameters are coupled (if not, they are read from file) 33 LOGICAL, PUBLIC :: cpl_hsig=.FALSE. 34 LOGICAL, PUBLIC :: cpl_phioc=.FALSE. 35 LOGICAL, PUBLIC :: cpl_sdrftx=.FALSE. 36 LOGICAL, PUBLIC :: cpl_sdrfty=.FALSE. 37 LOGICAL, PUBLIC :: cpl_wper=.FALSE. 38 LOGICAL, PUBLIC :: cpl_wnum=.FALSE. 39 LOGICAL, PUBLIC :: cpl_wstrf=.FALSE. 40 LOGICAL, PUBLIC :: cpl_wdrag=.FALSE. 41 42 INTEGER :: jpfld ! number of files to read for stokes drift 43 INTEGER :: jp_usd ! index of stokes drift (i-component) (m/s) at T-point 44 INTEGER :: jp_vsd ! index of stokes drift (j-component) (m/s) at T-point 45 INTEGER :: jp_swh ! index of significant wave hight (m) at T-point 46 INTEGER :: jp_wmp ! index of mean wave period (s) at T-point 35 47 36 48 TYPE(FLD), ALLOCATABLE, DIMENSION(:) :: sf_cd ! structure of input fields (file informations, fields read) Drag Coefficient … … 42 54 REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: tauoc_wave 43 55 REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: tsd2d 56 REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: usd2d, vsd2d 57 REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: zusd2dt, zvsd2dt 44 58 REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:,:) :: usd3d, vsd3d, wsd3d 45 REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:,:) :: usd3dt , vsd3dt59 REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:,:) :: usd3dt, vsd3dt 46 60 47 61 !! * Substitutions … … 54 68 !!---------------------------------------------------------------------- 55 69 CONTAINS 70 71 SUBROUTINE sbc_stokes( ) 72 !!--------------------------------------------------------------------- 73 !! *** ROUTINE sbc_stokes *** 74 !! 75 !! ** Purpose : compute the 3d Stokes Drift according to Breivik et al., 76 !! 2014 (DOI: 10.1175/JPO-D-14-0020.1) 77 !! 78 !! ** Method : - Calculate Stokes transport speed 79 !! - Calculate horizontal divergence 80 !! - Integrate the horizontal divergenze from the bottom 81 !! ** action 82 !!--------------------------------------------------------------------- 83 INTEGER :: jj,ji,jk 84 REAL(wp) :: ztransp, zfac, zsp0, zk, zus, zvs 85 REAL(wp), DIMENSION(:,:,:), POINTER :: ze3hdiv ! 3D workspace 86 !!--------------------------------------------------------------------- 87 ! 88 89 CALL wrk_alloc( jpi,jpj,jpk, ze3hdiv ) 90 DO jk = 1, jpk 91 DO jj = 1, jpj 92 DO ji = 1, jpi 93 ! On T grid 94 ! Stokes transport speed estimated from Hs and Tmean 95 ztransp = 2.0_wp*rpi*swh(ji,jj)**2.0_wp/(16.0_wp*MAX(wmp(ji,jj),0.0000001_wp)) 96 ! Stokes surface speed 97 zsp0 = SQRT( zusd2dt(ji,jj)**2 + zvsd2dt(ji,jj)**2) 98 ! Wavenumber scale 99 zk = ABS(zsp0)/MAX(ABS(5.97_wp*ztransp),0.0000001_wp) 100 ! Depth attenuation 101 zfac = EXP(-2.0_wp*zk*fsdept(ji,jj,jk))/(1.0_wp+8.0_wp*zk*fsdept(ji,jj,jk)) 102 ! 103 usd3dt(ji,jj,jk) = zfac * zusd2dt(ji,jj) * tmask(ji,jj,jk) 104 vsd3dt(ji,jj,jk) = zfac * zvsd2dt(ji,jj) * tmask(ji,jj,jk) 105 END DO 106 END DO 107 END DO 108 ! Into the U and V Grid 109 DO jk = 1, jpkm1 110 DO jj = 1, jpjm1 111 DO ji = 1, fs_jpim1 112 usd3d(ji,jj,jk) = 0.5 * umask(ji,jj,jk) * & 113 & ( usd3dt(ji,jj,jk) + usd3dt(ji+1,jj,jk) ) 114 vsd3d(ji,jj,jk) = 0.5 * vmask(ji,jj,jk) * & 115 & ( vsd3dt(ji,jj,jk) + vsd3dt(ji,jj+1,jk) ) 116 END DO 117 END DO 118 END DO 119 ! 120 CALL lbc_lnk( usd3d(:,:,:), 'U', -1. ) 121 CALL lbc_lnk( vsd3d(:,:,:), 'V', -1. ) 122 ! 123 DO jk = 1, jpkm1 ! Horizontal divergence 124 DO jj = 2, jpj 125 DO ji = fs_2, jpi 126 ze3hdiv(ji,jj,jk) = ( e2u(ji ,jj) * usd3d(ji ,jj,jk) & 127 & - e2u(ji-1,jj) * usd3d(ji-1,jj,jk) & 128 & + e1v(ji,jj ) * vsd3d(ji,jj ,jk) & 129 & - e1v(ji,jj-1) * vsd3d(ji,jj-1,jk) ) * r1_e1e2t(ji,jj) 130 END DO 131 END DO 132 END DO 133 ! 134 IF( .NOT. AGRIF_Root() ) THEN 135 IF( nbondi == 1 .OR. nbondi == 2 ) ze3hdiv(nlci-1, : ,:) = 0._wp ! east 136 IF( nbondi == -1 .OR. nbondi == 2 ) ze3hdiv( 2 , : ,:) = 0._wp ! west 137 IF( nbondj == 1 .OR. nbondj == 2 ) ze3hdiv( : ,nlcj-1,:) = 0._wp ! north 138 IF( nbondj == -1 .OR. nbondj == 2 ) ze3hdiv( : , 2 ,:) = 0._wp ! south 139 ENDIF 140 ! 141 CALL lbc_lnk( ze3hdiv, 'T', 1. ) 142 ! 143 DO jk = jpkm1, 1, -1 ! integrate from the bottom the e3t * hor. divergence 144 wsd3d(:,:,jk) = wsd3d(:,:,jk+1) - fse3t_n(:,:,jk) * ze3hdiv(:,:,jk) 145 END DO 146 #if defined key_bdy 147 IF( lk_bdy ) THEN 148 DO jk = 1, jpkm1 149 wsd3d(:,:,jk) = wsd3d(:,:,jk) * bdytmask(:,:) 150 END DO 151 ENDIF 152 #endif 153 CALL wrk_dealloc( jpi,jpj,jpk, ze3hdiv ) 154 ! 155 END SUBROUTINE sbc_stokes 156 157 SUBROUTINE sbc_qiao 158 !!--------------------------------------------------------------------- 159 !! *** ROUTINE sbc_qiao *** 160 !! 161 !! ** Purpose : Qiao formulation for wave enhanced turbulence 162 !! 2010 (DOI: 10.1007/s10236-010-0326) 163 !! 164 !! ** Method : - 165 !! ** action 166 !!--------------------------------------------------------------------- 167 INTEGER :: jj, ji 168 169 ! Calculate the module of the stokes drift on T grid 170 !------------------------------------------------- 171 DO jj = 1, jpj 172 DO ji = 1, jpi 173 tsd2d(ji,jj) = SQRT( zusd2dt(ji,jj) * zusd2dt(ji,jj) + zvsd2dt(ji,jj) * zvsd2dt(ji,jj) ) 174 END DO 175 END DO 176 ! 177 END SUBROUTINE sbc_qiao 56 178 57 179 SUBROUTINE sbc_wave( kt ) … … 73 195 INTEGER, INTENT( in ) :: kt ! ocean time step 74 196 ! 75 INTEGER :: ierror ! return error code 76 INTEGER :: ifpr, jj,ji,jk ! dummy loop indice 77 INTEGER :: ios ! Local integer output status for namelist read 78 CHARACTER(len=100) :: cn_dir ! Root directory for location of drag coefficient files 79 REAL(wp) :: ztransp, zfac, zsp0, zk, zus, zvs 80 REAL(wp), DIMENSION(jpi,jpj) :: zusd2dt, zvsd2dt 81 REAL(wp), DIMENSION(jpi,jpj,jpk) :: ze3hdiv ! 3D workspace 82 TYPE(FLD_N) :: sn_cdg, sn_usd, sn_vsd ! informations about the fields to be read 83 TYPE(FLD_N) :: sn_swh, sn_wmp, sn_wnum, sn_tauoc ! " " " " " " " 84 TYPE(FLD_N), DIMENSION(jpfld) :: slf_i ! array of namelist informations on the fields to read 197 INTEGER :: ierror ! return error code 198 INTEGER :: ifpr 199 INTEGER :: ios ! Local integer output status for namelist read 200 ! 201 CHARACTER(len=100) :: cn_dir ! Root directory for location of drag coefficient files 202 TYPE(FLD_N), ALLOCATABLE, DIMENSION(:) :: slf_i ! array of namelist informations on the fields to read 203 TYPE(FLD_N) :: sn_cdg, sn_usd, sn_vsd, & 204 & sn_swh, sn_wmp, sn_wnum, sn_tauoc ! informations about the fields to be read 205 !! 85 206 NAMELIST/namsbc_wave/ sn_cdg, cn_dir, sn_usd, sn_vsd, sn_swh, sn_wmp, sn_wnum, sn_tauoc 86 207 !!--------------------------------------------------------------------- … … 99 220 ! 100 221 IF( ln_cdgw ) THEN 101 ALLOCATE( sf_cd(1), STAT=ierror ) ! allocate and fill sf_wave with sn_cdg 102 IF( ierror > 0 ) CALL ctl_stop( 'STOP', 'sbc_wave: unable to allocate sf_wave structure' ) 103 ! 104 ALLOCATE( sf_cd(1)%fnow(jpi,jpj,1) ) 105 IF( sn_cdg%ln_tint ) ALLOCATE( sf_cd(1)%fdta(jpi,jpj,1,2) ) 106 CALL fld_fill( sf_cd, (/ sn_cdg /), cn_dir, 'sbc_wave', 'Wave module ', 'namsbc_wave' ) 222 IF( .NOT. cpl_wdrag ) THEN 223 ALLOCATE( sf_cd(1), STAT=ierror ) !* allocate and fill sf_wave with sn_cdg 224 IF( ierror > 0 ) CALL ctl_stop( 'STOP', 'sbc_wave: unable to allocate sf_wave structure' ) 225 ! 226 ALLOCATE( sf_cd(1)%fnow(jpi,jpj,1) ) 227 IF( sn_cdg%ln_tint ) ALLOCATE( sf_cd(1)%fdta(jpi,jpj,1,2) ) 228 CALL fld_fill( sf_cd, (/ sn_cdg /), cn_dir, 'sbc_wave', 'Wave module ', 'namsbc_wave' ) 229 ENDIF 107 230 ALLOCATE( cdn_wave(jpi,jpj) ) 108 231 ENDIF 109 232 110 233 IF( ln_tauoc ) THEN 111 ALLOCATE( sf_tauoc(1), STAT=ierror ) ! allocate and fill sf_wave with sn_tauoc 112 IF( ierror > 0 ) CALL ctl_stop( 'STOP', 'sbc_wave: unable to allocate sf_wave structure' ) 113 ! 114 ALLOCATE( sf_tauoc(1)%fnow(jpi,jpj,1) ) 115 IF( sn_tauoc%ln_tint ) ALLOCATE( sf_tauoc(1)%fdta(jpi,jpj,1,2) ) 116 CALL fld_fill( sf_tauoc, (/ sn_tauoc /), cn_dir, 'sbc_wave', 'Wave module', 'namsbc_wave' ) 234 IF( .NOT. cpl_wstrf ) THEN 235 ALLOCATE( sf_tauoc(1), STAT=ierror ) !* allocate and fill sf_wave with sn_tauoc 236 IF( ierror > 0 ) CALL ctl_stop( 'STOP', 'sbc_wave: unable to allocate sf_wave structure' ) 237 ! 238 ALLOCATE( sf_tauoc(1)%fnow(jpi,jpj,1) ) 239 IF( sn_tauoc%ln_tint ) ALLOCATE( sf_tauoc(1)%fdta(jpi,jpj,1,2) ) 240 CALL fld_fill( sf_tauoc, (/ sn_tauoc /), cn_dir, 'sbc_wave', 'Wave module', 'namsbc_wave' ) 241 ENDIF 117 242 ALLOCATE( tauoc_wave(jpi,jpj) ) 118 tauoc_wave(:,:) = 0.0119 243 ENDIF 120 244 121 245 IF( ln_sdw ) THEN 122 slf_i(jp_usd) = sn_usd ; slf_i(jp_vsd) = sn_vsd; 123 slf_i(jp_swh) = sn_swh ; slf_i(jp_wmp) = sn_wmp; 124 ALLOCATE( sf_sd(jpfld), STAT=ierror ) ! allocate and fill sf_sd with stokes drift 125 IF( ierror > 0 ) CALL ctl_stop( 'STOP', 'sbc_wave: unable to allocate sf_wave structure' ) 126 ! 127 DO ifpr = 1, jpfld 128 ALLOCATE( sf_sd(ifpr)%fnow(jpi,jpj,1) ) 129 IF( slf_i(ifpr)%ln_tint ) ALLOCATE( sf_sd(ifpr)%fdta(jpi,jpj,1,2) ) 130 END DO 131 132 CALL fld_fill( sf_sd, slf_i, cn_dir, 'sbc_wave', 'Wave module ', 'namsbc_wave' ) 133 ALLOCATE( usd3dt(jpi,jpj,jpk), vsd3dt(jpi,jpj,jpk), wsd3d(jpi,jpj,jpk) ) 134 ALLOCATE( usd3d (jpi,jpj,jpk), vsd3d (jpi,jpj,jpk) ) 246 ! Find out how many fields have to be read from file if not coupled 247 jpfld=0 248 jp_usd=0; jp_vsd=0; jp_swh=0; jp_wmp=0 249 IF( .NOT. cpl_sdrftx ) THEN 250 jpfld=jpfld+1 251 jp_usd=jpfld 252 ENDIF 253 IF( .NOT. cpl_sdrfty ) THEN 254 jpfld=jpfld+1 255 jp_vsd=jpfld 256 ENDIF 257 IF( .NOT. cpl_hsig ) THEN 258 jpfld=jpfld+1 259 jp_swh=jpfld 260 ENDIF 261 IF( .NOT. cpl_wper ) THEN 262 jpfld=jpfld+1 263 jp_wmp=jpfld 264 ENDIF 265 266 ! Read from file only the non-coupled fields 267 IF( jpfld > 0 ) THEN 268 ALLOCATE( slf_i(jpfld) ) 269 IF( jp_usd > 0 ) slf_i(jp_usd) = sn_usd 270 IF( jp_vsd > 0 ) slf_i(jp_vsd) = sn_vsd 271 IF( jp_swh > 0 ) slf_i(jp_swh) = sn_swh 272 IF( jp_wmp > 0 ) slf_i(jp_wmp) = sn_wmp 273 ALLOCATE( sf_sd(jpfld), STAT=ierror ) !* allocate and fill sf_sd with stokes drift 274 IF( ierror > 0 ) CALL ctl_stop( 'STOP', 'sbc_wave: unable to allocate sf_wave structure' ) 275 ! 276 DO ifpr= 1, jpfld 277 ALLOCATE( sf_sd(ifpr)%fnow(jpi,jpj,1) ) 278 IF( slf_i(ifpr)%ln_tint ) ALLOCATE( sf_sd(ifpr)%fdta(jpi,jpj,1,2) ) 279 END DO 280 281 CALL fld_fill( sf_sd, slf_i, cn_dir, 'sbc_wave', 'Wave module ', 'namsbc_wave' ) 282 ENDIF 283 ALLOCATE( usd2d(jpi,jpj),vsd2d(jpi,jpj) ) 284 ALLOCATE( usd3d(jpi,jpj,jpk),vsd3d(jpi,jpj,jpk),wsd3d(jpi,jpj,jpk) ) 285 ALLOCATE( usd3dt(jpi,jpj,jpk),vsd3dt(jpi,jpj,jpk) ) 135 286 ALLOCATE( swh(jpi,jpj), wmp(jpi,jpj) ) 136 usd3d(:,:,:) = 0._wp 137 vsd3d(:,:,:) = 0._wp 138 wsd3d(:,:,:) = 0._wp 287 ALLOCATE( zusd2dt(jpi,jpj), zvsd2dt(jpi,jpj) ) 288 usd3d(:,:,:) = 0._wp 289 vsd3d(:,:,:) = 0._wp 290 wsd3d(:,:,:) = 0._wp 139 291 IF( ln_zdfqiao ) THEN !== Vertical mixing enhancement using Qiao,2010 ==! 140 ALLOCATE( sf_wn(1), STAT=ierror ) !* allocate and fill sf_wave with sn_wnum 141 IF( ierror > 0 ) CALL ctl_stop( 'STOP', 'sbc_wave: unable to allocate sf_wave structure' ) 142 ALLOCATE( sf_wn(1)%fnow(jpi,jpj,1) ) 143 IF( sn_wnum%ln_tint ) ALLOCATE( sf_wn(1)%fdta(jpi,jpj,1,2) ) 144 CALL fld_fill( sf_wn, (/ sn_wnum /), cn_dir, 'sbc_wave', 'Wave module', 'namsbc_wave' ) 292 IF( .NOT. cpl_wnum ) THEN 293 ALLOCATE( sf_wn(1), STAT=ierror ) !* allocate and fill sf_wave with sn_wnum 294 IF( ierror > 0 ) CALL ctl_stop( 'STOP', 'sbc_wave: unable toallocate sf_wave structure' ) 295 ALLOCATE( sf_wn(1)%fnow(jpi,jpj,1) ) 296 IF( sn_wnum%ln_tint ) ALLOCATE( sf_wn(1)%fdta(jpi,jpj,1,2) ) 297 CALL fld_fill( sf_wn, (/ sn_wnum /), cn_dir, 'sbc_wave', 'Wave module', 'namsbc_wave' ) 298 ENDIF 145 299 ALLOCATE( wnum(jpi,jpj),tsd2d(jpi,jpj) ) 146 300 ENDIF … … 148 302 ENDIF 149 303 ! 150 IF( ln_cdgw ) THEN !== Neutral drag coefficient ==!151 CALL fld_read( kt, nn_fsbc, sf_cd ) 304 IF( ln_cdgw .AND. .NOT. cpl_wdrag ) THEN !== Neutral drag coefficient ==! 305 CALL fld_read( kt, nn_fsbc, sf_cd ) ! read from external forcing 152 306 cdn_wave(:,:) = sf_cd(1)%fnow(:,:,1) 153 307 ENDIF 154 308 155 IF( ln_tauoc ) THEN !== Wave induced stress ==!156 CALL fld_read( kt, nn_fsbc, sf_tauoc ) !read wave norm stress from external forcing309 IF( ln_tauoc .AND. .NOT. cpl_wstrf ) THEN !== Wave induced stress ==! 310 CALL fld_read( kt, nn_fsbc, sf_tauoc ) !* read wave norm stress from external forcing 157 311 tauoc_wave(:,:) = sf_tauoc(1)%fnow(:,:,1) 158 312 ENDIF 159 313 160 IF( ln_sdw ) THEN !== Computation of the 3d Stokes Drift ==! 161 CALL fld_read( kt, nn_fsbc, sf_sd ) ! read wave parameters from external forcing 162 swh(:,:) = sf_sd(jp_swh)%fnow(:,:,1) ! significant wave height 163 wmp(:,:) = sf_sd(jp_wmp)%fnow(:,:,1) ! wave mean period 164 zusd2dt(:,:) = sf_sd(jp_usd)%fnow(:,:,1) ! 2D zonal Stokes Drift at T point 165 zvsd2dt(:,:) = sf_sd(jp_vsd)%fnow(:,:,1) ! 2D meridional Stokes Drift at T point 314 IF( ln_sdw ) THEN !== Computation of the 3d Stokes Drift ==! 166 315 ! 316 ! Read from file only if the field is not coupled 317 IF( jpfld > 0 ) THEN 318 CALL fld_read( kt, nn_fsbc, sf_sd ) !* read wave parameters from external forcing 319 IF( jp_swh > 0 ) swh(:,:) = sf_sd(jp_swh)%fnow(:,:,1) ! significant wave height 320 IF( jp_wmp > 0 ) wmp(:,:) = sf_sd(jp_wmp)%fnow(:,:,1) ! wave mean period 321 IF( jp_usd > 0 ) zusd2dt(:,:) = sf_sd(jp_usd)%fnow(:,:,1) ! 2D zonal Stokes Drift at T point 322 IF( jp_vsd > 0 ) zvsd2dt(:,:) = sf_sd(jp_vsd)%fnow(:,:,1) ! 2D meridional Stokes Drift at T point 323 ENDIF 324 ! 325 ! Read also wave number if needed, so that it is available in coupling routines 326 IF( ln_zdfqiao .AND. .NOT. cpl_wnum ) THEN 327 CALL fld_read( kt, nn_fsbc, sf_wn ) !* read wave parameters from external forcing 328 wnum(:,:) = sf_wn(1)%fnow(:,:,1) 329 ENDIF 330 167 331 !== Computation of the 3d Stokes Drift according to Breivik et al.,2014 168 332 !(DOI: 10.1175/JPO-D-14-0020.1)==! 169 333 ! 170 DO jk = 1, jpk 171 DO jj = 1, jpj 172 DO ji = 1, jpi 173 ! On T grid 174 ! Stokes transport speed estimated from Hs and Tmean 175 ztransp = 2.0_wp*rpi*swh(ji,jj)**2.0_wp/(16.0_wp*MAX(wmp(ji,jj),0.0000001_wp)) 176 ! Stokes surface speed 177 zsp0 = SQRT( zusd2dt(ji,jj)**2 + zvsd2dt(ji,jj)**2 ) 178 ! Wavenumber scale 179 zk = ABS(zsp0)/MAX(ABS(5.97_wp*ztransp),0.0000001_wp) 180 ! Depth attenuation 181 zfac = EXP(-2.0_wp*zk*fsdept(ji,jj,jk))/(1.0_wp+8.0_wp*zk*fsdept(ji,jj,jk)) 182 ! 183 usd3dt(ji,jj,jk) = zfac * zusd2dt(ji,jj) * tmask(ji,jj,jk) 184 vsd3dt(ji,jj,jk) = zfac * zvsd2dt(ji,jj) * tmask(ji,jj,jk) 185 END DO 186 END DO 187 END DO 188 ! Into the U and V Grid 189 DO jk = 1, jpkm1 190 DO jj = 1, jpjm1 191 DO ji = 1, jpim1 192 usd3d(ji,jj,jk) = 0.5 * umask(ji,jj,jk) * & 193 & ( usd3dt(ji,jj,jk) + usd3dt(ji+1,jj,jk) ) 194 vsd3d(ji,jj,jk) = 0.5 * vmask(ji,jj,jk) * & 195 & ( vsd3dt(ji,jj,jk) + vsd3dt(ji,jj+1,jk) ) 196 END DO 197 END DO 198 END DO 199 ! 200 CALL lbc_lnk( usd3d(:,:,:), 'U', -1. ) 201 CALL lbc_lnk( vsd3d(:,:,:), 'V', -1. ) 202 ! 203 DO jk = 1, jpkm1 ! Horizontal divergence 204 DO jj = 2, jpj 205 DO ji = fs_2, jpi 206 ze3hdiv(ji,jj,jk) = ( e2u(ji ,jj) * usd3d(ji ,jj,jk) & 207 & - e2u(ji-1,jj) * usd3d(ji-1,jj,jk) & 208 & + e1v(ji,jj ) * vsd3d(ji,jj ,jk) & 209 & - e1v(ji,jj-1) * vsd3d(ji,jj-1,jk) ) * r1_e1e2t(ji,jj) 210 END DO 211 END DO 212 END DO 213 ! 214 IF( .NOT. AGRIF_Root() ) THEN 215 IF( nbondi == 1 .OR. nbondi == 2 ) ze3hdiv(nlci-1, : ,:) = 0._wp ! east 216 IF( nbondi == -1 .OR. nbondi == 2 ) ze3hdiv( 2 , : ,:) = 0._wp ! west 217 IF( nbondj == 1 .OR. nbondj == 2 ) ze3hdiv( : ,nlcj-1,:) = 0._wp ! north 218 IF( nbondj == -1 .OR. nbondj == 2 ) ze3hdiv( : , 2 ,:) = 0._wp ! south 219 ENDIF 220 ! 221 CALL lbc_lnk( ze3hdiv, 'T', 1. ) 222 ! 223 DO jk = jpkm1, 1, -1 ! integrate from the bottom the e3t * hor. divergence 224 wsd3d(:,:,jk) = wsd3d(:,:,jk+1) - fse3t_n(:,:,jk) * ze3hdiv(:,:,jk) 225 END DO 226 #if defined key_bdy 227 IF( lk_bdy ) THEN 228 DO jk = 1, jpkm1 229 wsd3d(:,:,jk) = wsd3d(:,:,jk) * bdytmask(:,:) 230 END DO 231 ENDIF 232 #endif 233 234 IF ( ln_zdfqiao ) THEN 235 CALL fld_read( kt, nn_fsbc, sf_wn ) ! read wave parameters from external forcing 236 wnum(:,:) = sf_wn(1)%fnow(:,:,1) 237 238 ! Calculate the module of the stokes drift on T grid 239 !------------------------------------------------- 240 DO jj = 1, jpj 241 DO ji = 1, jpi 242 tsd2d(ji,jj) = SQRT( zusd2dt(ji,jj)**2 + zvsd2dt(ji,jj)**2 ) 243 END DO 244 END DO 245 ENDIF 246 ! 334 ! Calculate only if no necessary fields are coupled, if not calculate later after coupling 335 IF( jpfld == 4 ) THEN 336 CALL sbc_stokes() 337 IF( ln_zdfqiao .AND. .NOT. cpl_wnum ) THEN 338 CALL sbc_qiao() 339 ENDIF 340 ENDIF 247 341 ENDIF 248 342 ! -
branches/2015/dev_r5936_INGV1_WAVE/NEMOGCM/NEMO/OPA_SRC/step.F90
r7340 r7343 220 220 CALL dyn_adv ( kstp ) ! advection (vector or flux form) 221 221 CALL dyn_vor ( kstp ) ! vorticity term including Coriolis 222 IF( ln_wave .AND. ln_sdw .AND. ln_stcor ) & 223 & CALL dyn_stcor ( kstp ) ! Stokes-Coriolis forcing 222 224 CALL dyn_ldf ( kstp ) ! lateral mixing 223 225 IF( ln_wave .AND. ln_sdw .AND. ln_stcor) &
Note: See TracChangeset
for help on using the changeset viewer.