New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
sbcwave.F90 in NEMO/branches/2020/dev_r13383_HPC-02_Daley_Tiling/src/OCE/SBC – NEMO

source: NEMO/branches/2020/dev_r13383_HPC-02_Daley_Tiling/src/OCE/SBC/sbcwave.F90 @ 13553

Last change on this file since 13553 was 13553, checked in by hadcv, 4 years ago

Merge in trunk up to [13550]

  • Property svn:keywords set to Id
File size: 24.2 KB
RevLine 
[2990]1MODULE sbcwave
2   !!======================================================================
3   !!                       ***  MODULE  sbcwave  ***
4   !! Wave module
5   !!======================================================================
[7646]6   !! History :  3.3  !  2011-09  (M. Adani)  Original code: Drag Coefficient
7   !!         :  3.4  !  2012-10  (M. Adani)  Stokes Drift
8   !!            3.6  !  2014-09  (E. Clementi,P. Oddo) New Stokes Drift Computation
9   !!             -   !  2016-12  (G. Madec, E. Clementi) update Stoke drift computation
10   !!                                                    + add sbc_wave_ini routine
[2990]11   !!----------------------------------------------------------------------
12
13   !!----------------------------------------------------------------------
[7646]14   !!   sbc_stokes    : calculate 3D Stokes-drift velocities
15   !!   sbc_wave      : wave data from wave model in netcdf files
16   !!   sbc_wave_init : initialisation fo surface waves
[2990]17   !!----------------------------------------------------------------------
[7646]18   USE phycst         ! physical constants
19   USE oce            ! ocean variables
[5836]20   USE sbc_oce        ! Surface boundary condition: ocean fields
[9019]21   USE zdf_oce,  ONLY : ln_zdfswm
[7646]22   USE bdy_oce        ! open boundary condition variables
23   USE domvvl         ! domain: variable volume layers
[5836]24   !
25   USE iom            ! I/O manager library
26   USE in_out_manager ! I/O manager
27   USE lib_mpp        ! distribued memory computing library
28   USE fldread        ! read input fields
[2990]29
30   IMPLICIT NONE
31   PRIVATE
32
[7646]33   PUBLIC   sbc_stokes      ! routine called in sbccpl
[9023]34   PUBLIC   sbc_wstress     ! routine called in sbcmod
[7646]35   PUBLIC   sbc_wave        ! routine called in sbcmod
36   PUBLIC   sbc_wave_init   ! routine called in sbcmod
[2990]37   
[7646]38   ! Variables checking if the wave parameters are coupled (if not, they are read from file)
39   LOGICAL, PUBLIC ::   cpl_hsig   = .FALSE.
40   LOGICAL, PUBLIC ::   cpl_phioc  = .FALSE.
41   LOGICAL, PUBLIC ::   cpl_sdrftx = .FALSE.
42   LOGICAL, PUBLIC ::   cpl_sdrfty = .FALSE.
43   LOGICAL, PUBLIC ::   cpl_wper   = .FALSE.
[9023]44   LOGICAL, PUBLIC ::   cpl_wfreq  = .FALSE.
[7646]45   LOGICAL, PUBLIC ::   cpl_wnum   = .FALSE.
[9115]46   LOGICAL, PUBLIC ::   cpl_tauwoc = .FALSE.
[9023]47   LOGICAL, PUBLIC ::   cpl_tauw   = .FALSE.
[7646]48   LOGICAL, PUBLIC ::   cpl_wdrag  = .FALSE.
[5836]49
[7646]50   INTEGER ::   jpfld    ! number of files to read for stokes drift
51   INTEGER ::   jp_usd   ! index of stokes drift  (i-component) (m/s)    at T-point
52   INTEGER ::   jp_vsd   ! index of stokes drift  (j-component) (m/s)    at T-point
53   INTEGER ::   jp_hsw   ! index of significant wave hight      (m)      at T-point
54   INTEGER ::   jp_wmp   ! index of mean wave period            (s)      at T-point
[9023]55   INTEGER ::   jp_wfr   ! index of wave peak frequency         (1/s)    at T-point
[2990]56
[7646]57   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf_cd      ! structure of input fields (file informations, fields read) Drag Coefficient
58   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf_sd      ! structure of input fields (file informations, fields read) Stokes Drift
59   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf_wn      ! structure of input fields (file informations, fields read) wave number for Qiao
[9115]60   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf_tauwoc  ! structure of input fields (file informations, fields read) normalized wave stress into the ocean
[9023]61   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf_tauw    ! structure of input fields (file informations, fields read) ocean stress components from wave model
62
[7646]63   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)   ::   cdn_wave            !:
64   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)   ::   hsw, wmp, wnum      !:
[9023]65   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)   ::   wfreq               !:
[7646]66   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)   ::   tauoc_wave          !: 
[9023]67   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)   ::   tauw_x, tauw_y      !: 
[7646]68   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)   ::   tsd2d               !:
69   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)   ::   div_sd              !: barotropic stokes drift divergence
70   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)   ::   ut0sd, vt0sd        !: surface Stokes drift velocities at t-point
71   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:,:) ::   usd  , vsd  , wsd   !: Stokes drift velocities at u-, v- & w-points, resp.
[5836]72
[3680]73   !! * Substitutions
[12377]74#  include "do_loop_substitute.h90"
[13237]75#  include "domzgr_substitute.h90"
[2990]76   !!----------------------------------------------------------------------
[10068]77   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
[5215]78   !! $Id$
[10068]79   !! Software governed by the CeCILL license (see ./LICENSE)
[2990]80   !!----------------------------------------------------------------------
81CONTAINS
82
[12377]83   SUBROUTINE sbc_stokes( Kmm )
[7646]84      !!---------------------------------------------------------------------
85      !!                     ***  ROUTINE sbc_stokes  ***
86      !!
87      !! ** Purpose :   compute the 3d Stokes Drift according to Breivik et al.,
88      !!                2014 (DOI: 10.1175/JPO-D-14-0020.1)
89      !!
90      !! ** Method  : - Calculate Stokes transport speed
91      !!              - Calculate horizontal divergence
92      !!              - Integrate the horizontal divergenze from the bottom
93      !! ** action 
94      !!---------------------------------------------------------------------
[12377]95      INTEGER, INTENT(in) :: Kmm ! ocean time level index
[7646]96      INTEGER  ::   jj, ji, jk   ! dummy loop argument
97      INTEGER  ::   ik           ! local integer
[9019]98      REAL(wp) ::  ztransp, zfac, zsp0
99      REAL(wp) ::  zdepth, zsqrt_depth,  zexp_depth, z_two_thirds, zsqrtpi !sqrt of pi
100      REAL(wp) ::  zbot_u, zbot_v, zkb_u, zkb_v, zke3_u, zke3_v, zda_u, zda_v
101      REAL(wp) ::  zstokes_psi_u_bot, zstokes_psi_v_bot
[9029]102      REAL(wp) ::  zdep_u, zdep_v, zkh_u, zkh_v
[9115]103      REAL(wp), DIMENSION(:,:)  , ALLOCATABLE ::   zk_t, zk_u, zk_v, zu0_sd, zv0_sd     ! 2D workspace
104      REAL(wp), DIMENSION(:,:)  , ALLOCATABLE ::   zstokes_psi_u_top, zstokes_psi_v_top ! 2D workspace
105      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::   ze3divh                              ! 3D workspace
[7646]106      !!---------------------------------------------------------------------
107      !
[13553]108      ALLOCATE( ze3divh(jpi,jpj,jpkm1) )   ! jpkm1 -> avoid lbc_lnk on jpk that is not defined
[9115]109      ALLOCATE( zk_t(jpi,jpj), zk_u(jpi,jpj), zk_v(jpi,jpj), zu0_sd(jpi,jpj), zv0_sd(jpi,jpj) )
[7646]110      !
[9023]111      ! select parameterization for the calculation of vertical Stokes drift
112      ! exp. wave number at t-point
[9115]113      IF( ll_st_bv_li ) THEN   ! (Eq. (19) in Breivik et al. (2014) )
[9023]114         zfac = 2.0_wp * rpi / 16.0_wp
[13295]115         DO_2D( 1, 1, 1, 1 )
[12377]116            ! Stokes drift velocity estimated from Hs and Tmean
117            ztransp = zfac * hsw(ji,jj)*hsw(ji,jj) / MAX( wmp(ji,jj), 0.0000001_wp )
118            ! Stokes surface speed
119            tsd2d(ji,jj) = SQRT( ut0sd(ji,jj)*ut0sd(ji,jj) + vt0sd(ji,jj)*vt0sd(ji,jj))
120            ! Wavenumber scale
121            zk_t(ji,jj) = ABS( tsd2d(ji,jj) ) / MAX( ABS( 5.97_wp*ztransp ), 0.0000001_wp )
122         END_2D
[13553]123         DO_2D( 1, 0, 1, 0 )          ! exp. wave number & Stokes drift velocity at u- & v-points
[12377]124            zk_u(ji,jj) = 0.5_wp * ( zk_t(ji,jj) + zk_t(ji+1,jj) )
125            zk_v(ji,jj) = 0.5_wp * ( zk_t(ji,jj) + zk_t(ji,jj+1) )
126            !
127            zu0_sd(ji,jj) = 0.5_wp * ( ut0sd(ji,jj) + ut0sd(ji+1,jj) )
128            zv0_sd(ji,jj) = 0.5_wp * ( vt0sd(ji,jj) + vt0sd(ji,jj+1) )
129         END_2D
[9115]130      ELSE IF( ll_st_peakfr ) THEN    ! peak wave number calculated from the peak frequency received by the wave model
[13295]131         DO_2D( 1, 1, 1, 1 )
[12377]132            zk_t(ji,jj) = ( 2.0_wp * rpi * wfreq(ji,jj) ) * ( 2.0_wp * rpi * wfreq(ji,jj) ) / grav
133         END_2D
[13295]134         DO_2D( 1, 0, 1, 0 )
[12377]135            zk_u(ji,jj) = 0.5_wp * ( zk_t(ji,jj) + zk_t(ji+1,jj) )
136            zk_v(ji,jj) = 0.5_wp * ( zk_t(ji,jj) + zk_t(ji,jj+1) )
137            !
138            zu0_sd(ji,jj) = 0.5_wp * ( ut0sd(ji,jj) + ut0sd(ji+1,jj) )
139            zv0_sd(ji,jj) = 0.5_wp * ( vt0sd(ji,jj) + vt0sd(ji,jj+1) )
140         END_2D
[9023]141      ENDIF
[7646]142      !
143      !                       !==  horizontal Stokes Drift 3D velocity  ==!
[9115]144      IF( ll_st_bv2014 ) THEN
[13295]145         DO_3D( 0, 0, 0, 0, 1, jpkm1 )
[12377]146            zdep_u = 0.5_wp * ( gdept(ji,jj,jk,Kmm) + gdept(ji+1,jj,jk,Kmm) )
147            zdep_v = 0.5_wp * ( gdept(ji,jj,jk,Kmm) + gdept(ji,jj+1,jk,Kmm) )
148            !                         
149            zkh_u = zk_u(ji,jj) * zdep_u     ! k * depth
150            zkh_v = zk_v(ji,jj) * zdep_v
151            !                                ! Depth attenuation
152            zda_u = EXP( -2.0_wp*zkh_u ) / ( 1.0_wp + 8.0_wp*zkh_u )
153            zda_v = EXP( -2.0_wp*zkh_v ) / ( 1.0_wp + 8.0_wp*zkh_v )
154            !
155            usd(ji,jj,jk) = zda_u * zu0_sd(ji,jj) * umask(ji,jj,jk)
156            vsd(ji,jj,jk) = zda_v * zv0_sd(ji,jj) * vmask(ji,jj,jk)
157         END_3D
[9115]158      ELSE IF( ll_st_li2017 .OR. ll_st_peakfr ) THEN
159         ALLOCATE( zstokes_psi_u_top(jpi,jpj), zstokes_psi_v_top(jpi,jpj) )
[13295]160         DO_2D( 1, 0, 1, 0 )
[12377]161            zstokes_psi_u_top(ji,jj) = 0._wp
162            zstokes_psi_v_top(ji,jj) = 0._wp
163         END_2D
[9117]164         zsqrtpi = SQRT(rpi)
165         z_two_thirds = 2.0_wp / 3.0_wp
[13553]166         DO_3D( 0, 0, 0, 0, 1, jpkm1 )       ! exp. wave number & Stokes drift velocity at u- & v-points
[12377]167            zbot_u = ( gdepw(ji,jj,jk+1,Kmm) + gdepw(ji+1,jj,jk+1,Kmm) )  ! 2 * bottom depth
168            zbot_v = ( gdepw(ji,jj,jk+1,Kmm) + gdepw(ji,jj+1,jk+1,Kmm) )  ! 2 * bottom depth
169            zkb_u  = zk_u(ji,jj) * zbot_u                             ! 2 * k * bottom depth
170            zkb_v  = zk_v(ji,jj) * zbot_v                             ! 2 * k * bottom depth
171            !
172            zke3_u = MAX(1.e-8_wp, 2.0_wp * zk_u(ji,jj) * e3u(ji,jj,jk,Kmm))     ! 2k * thickness
173            zke3_v = MAX(1.e-8_wp, 2.0_wp * zk_v(ji,jj) * e3v(ji,jj,jk,Kmm))     ! 2k * thickness
[9115]174
[12377]175            ! Depth attenuation .... do u component first..
176            zdepth      = zkb_u
177            zsqrt_depth = SQRT(zdepth)
178            zexp_depth  = EXP(-zdepth)
179            zstokes_psi_u_bot = 1.0_wp - zexp_depth  &
180                 &              - z_two_thirds * ( zsqrtpi*zsqrt_depth*zdepth*ERFC(zsqrt_depth) &
181                 &              + 1.0_wp - (1.0_wp + zdepth)*zexp_depth )
182            zda_u                    = ( zstokes_psi_u_bot - zstokes_psi_u_top(ji,jj) ) / zke3_u
183            zstokes_psi_u_top(ji,jj) =   zstokes_psi_u_bot
[9115]184
[12377]185            !         ... and then v component
186            zdepth      =zkb_v
187            zsqrt_depth = SQRT(zdepth)
188            zexp_depth  = EXP(-zdepth)
189            zstokes_psi_v_bot = 1.0_wp - zexp_depth  &
190                 &              - z_two_thirds * ( zsqrtpi*zsqrt_depth*zdepth*ERFC(zsqrt_depth) &
191                 &              + 1.0_wp - (1.0_wp + zdepth)*zexp_depth )
192            zda_v                    = ( zstokes_psi_v_bot - zstokes_psi_v_top(ji,jj) ) / zke3_v
193            zstokes_psi_v_top(ji,jj) =   zstokes_psi_v_bot
194            !
195            usd(ji,jj,jk) = zda_u * zu0_sd(ji,jj) * umask(ji,jj,jk)
196            vsd(ji,jj,jk) = zda_v * zv0_sd(ji,jj) * vmask(ji,jj,jk)
197         END_3D
[9115]198         DEALLOCATE( zstokes_psi_u_top, zstokes_psi_v_top )
[9023]199      ENDIF
200
[13226]201      CALL lbc_lnk_multi( 'sbcwave', usd, 'U', -1.0_wp, vsd, 'V', -1.0_wp )
[9019]202
[7646]203      !
204      !                       !==  vertical Stokes Drift 3D velocity  ==!
205      !
[13553]206      DO_3D( 0, 1, 0, 1, 1, jpkm1 )    ! Horizontal e3*divergence
[12377]207         ze3divh(ji,jj,jk) = (  e2u(ji  ,jj) * e3u(ji  ,jj,jk,Kmm) * usd(ji  ,jj,jk)    &
208            &                 - e2u(ji-1,jj) * e3u(ji-1,jj,jk,Kmm) * usd(ji-1,jj,jk)    &
209            &                 + e1v(ji,jj  ) * e3v(ji,jj  ,jk,Kmm) * vsd(ji,jj  ,jk)    &
[13237]210            &                 - e1v(ji,jj-1) * e3v(ji,jj-1,jk,Kmm) * vsd(ji,jj-1,jk)  ) &
211            &                * r1_e1e2t(ji,jj)
[12377]212      END_3D
[7646]213      !
[13226]214      CALL lbc_lnk( 'sbcwave', ze3divh, 'T', 1.0_wp )
[7646]215      !
216      IF( ln_linssh ) THEN   ;   ik = 1   ! none zero velocity through the sea surface
217      ELSE                   ;   ik = 2   ! w=0 at the surface (set one for all in sbc_wave_init)
218      ENDIF
219      DO jk = jpkm1, ik, -1          ! integrate from the bottom the hor. divergence (NB: at k=jpk w is always zero)
220         wsd(:,:,jk) = wsd(:,:,jk+1) - ze3divh(:,:,jk)
221      END DO
222      !
223      IF( ln_bdy ) THEN
224         DO jk = 1, jpkm1
225            wsd(:,:,jk) = wsd(:,:,jk) * bdytmask(:,:)
226         END DO
227      ENDIF
228      !                       !==  Horizontal divergence of barotropic Stokes transport  ==!
229      div_sd(:,:) = 0._wp
230      DO jk = 1, jpkm1                                 !
231        div_sd(:,:) = div_sd(:,:) + ze3divh(:,:,jk)
232      END DO
233      !
234      CALL iom_put( "ustokes",  usd  )
235      CALL iom_put( "vstokes",  vsd  )
236      CALL iom_put( "wstokes",  wsd  )
237      !
[9115]238      DEALLOCATE( ze3divh )
239      DEALLOCATE( zk_t, zk_u, zk_v, zu0_sd, zv0_sd )
[7646]240      !
241   END SUBROUTINE sbc_stokes
242
243
[9023]244   SUBROUTINE sbc_wstress( )
245      !!---------------------------------------------------------------------
246      !!                     ***  ROUTINE sbc_wstress  ***
247      !!
248      !! ** Purpose :   Updates the ocean momentum modified by waves
249      !!
250      !! ** Method  : - Calculate u,v components of stress depending on stress
251      !!                model
252      !!              - Calculate the stress module
253      !!              - The wind module is not modified by waves
254      !! ** action 
255      !!---------------------------------------------------------------------
256      INTEGER  ::   jj, ji   ! dummy loop argument
257      !
[9033]258      IF( ln_tauwoc ) THEN
[9023]259         utau(:,:) = utau(:,:)*tauoc_wave(:,:)
260         vtau(:,:) = vtau(:,:)*tauoc_wave(:,:)
261         taum(:,:) = taum(:,:)*tauoc_wave(:,:)
262      ENDIF
263      !
264      IF( ln_tauw ) THEN
[13295]265         DO_2D( 1, 0, 1, 0 )
[12377]266            ! Stress components at u- & v-points
267            utau(ji,jj) = 0.5_wp * ( tauw_x(ji,jj) + tauw_x(ji+1,jj) )
268            vtau(ji,jj) = 0.5_wp * ( tauw_y(ji,jj) + tauw_y(ji,jj+1) )
269            !
270            ! Stress module at t points
271            taum(ji,jj) = SQRT( tauw_x(ji,jj)*tauw_x(ji,jj) + tauw_y(ji,jj)*tauw_y(ji,jj) )
272         END_2D
[13226]273         CALL lbc_lnk_multi( 'sbcwave', utau(:,:), 'U', -1.0_wp , vtau(:,:), 'V', -1.0_wp , taum(:,:) , 'T', -1.0_wp )
[9023]274      ENDIF
275      !
276   END SUBROUTINE sbc_wstress
277
278
[12377]279   SUBROUTINE sbc_wave( kt, Kmm )
[2990]280      !!---------------------------------------------------------------------
[7646]281      !!                     ***  ROUTINE sbc_wave  ***
[2990]282      !!
[7646]283      !! ** Purpose :   read wave parameters from wave model  in netcdf files.
[2990]284      !!
285      !! ** Method  : - Read namelist namsbc_wave
286      !!              - Read Cd_n10 fields in netcdf files
[3680]287      !!              - Read stokes drift 2d in netcdf files
[7646]288      !!              - Read wave number in netcdf files
289      !!              - Compute 3d stokes drift using Breivik et al.,2014
290      !!                formulation
291      !! ** action 
[2990]292      !!---------------------------------------------------------------------
[7646]293      INTEGER, INTENT(in   ) ::   kt   ! ocean time step
[12377]294      INTEGER, INTENT(in   ) ::   Kmm  ! ocean time index
[2990]295      !!---------------------------------------------------------------------
296      !
[7646]297      IF( ln_cdgw .AND. .NOT. cpl_wdrag ) THEN     !==  Neutral drag coefficient  ==!
298         CALL fld_read( kt, nn_fsbc, sf_cd )             ! read from external forcing
[9821]299         cdn_wave(:,:) = sf_cd(1)%fnow(:,:,1) * tmask(:,:,1)
[7646]300      ENDIF
301
[9115]302      IF( ln_tauwoc .AND. .NOT. cpl_tauwoc ) THEN  !==  Wave induced stress  ==!
303         CALL fld_read( kt, nn_fsbc, sf_tauwoc )         ! read wave norm stress from external forcing
[9821]304         tauoc_wave(:,:) = sf_tauwoc(1)%fnow(:,:,1) * tmask(:,:,1)
[7646]305      ENDIF
306
[9023]307      IF( ln_tauw .AND. .NOT. cpl_tauw ) THEN      !==  Wave induced stress  ==!
308         CALL fld_read( kt, nn_fsbc, sf_tauw )           ! read ocean stress components from external forcing (T grid)
[9821]309         tauw_x(:,:) = sf_tauw(1)%fnow(:,:,1) * tmask(:,:,1)
310         tauw_y(:,:) = sf_tauw(2)%fnow(:,:,1) * tmask(:,:,1)
[9023]311      ENDIF
312
[7646]313      IF( ln_sdw )  THEN                           !==  Computation of the 3d Stokes Drift  ==!
[6140]314         !
[7646]315         IF( jpfld > 0 ) THEN                            ! Read from file only if the field is not coupled
316            CALL fld_read( kt, nn_fsbc, sf_sd )          ! read wave parameters from external forcing
[9821]317            IF( jp_hsw > 0 )   hsw  (:,:) = sf_sd(jp_hsw)%fnow(:,:,1) * tmask(:,:,1)  ! significant wave height
318            IF( jp_wmp > 0 )   wmp  (:,:) = sf_sd(jp_wmp)%fnow(:,:,1) * tmask(:,:,1)  ! wave mean period
319            IF( jp_wfr > 0 )   wfreq(:,:) = sf_sd(jp_wfr)%fnow(:,:,1) * tmask(:,:,1)  ! Peak wave frequency
320            IF( jp_usd > 0 )   ut0sd(:,:) = sf_sd(jp_usd)%fnow(:,:,1) * tmask(:,:,1)  ! 2D zonal Stokes Drift at T point
321            IF( jp_vsd > 0 )   vt0sd(:,:) = sf_sd(jp_vsd)%fnow(:,:,1) * tmask(:,:,1)  ! 2D meridional Stokes Drift at T point
[7646]322         ENDIF
[2990]323         !
[7646]324         ! Read also wave number if needed, so that it is available in coupling routines
[9019]325         IF( ln_zdfswm .AND. .NOT.cpl_wnum ) THEN
[7646]326            CALL fld_read( kt, nn_fsbc, sf_wn )          ! read wave parameters from external forcing
[9821]327            wnum(:,:) = sf_wn(1)%fnow(:,:,1) * tmask(:,:,1)
[6140]328         ENDIF
[7646]329           
[9115]330         ! Calculate only if required fields have been read
331         ! In coupled wave model-NEMO case the call is done after coupling
[6140]332         !
[9115]333         IF( ( ll_st_bv_li   .AND. jp_hsw>0 .AND. jp_wmp>0 .AND. jp_usd>0 .AND. jp_vsd>0 ) .OR. &
[12377]334           & ( ll_st_peakfr  .AND. jp_wfr>0 .AND. jp_usd>0 .AND. jp_vsd>0                ) ) CALL sbc_stokes( Kmm )
[6140]335         !
[7646]336      ENDIF
337      !
338   END SUBROUTINE sbc_wave
339
340
341   SUBROUTINE sbc_wave_init
342      !!---------------------------------------------------------------------
343      !!                     ***  ROUTINE sbc_wave_init  ***
344      !!
345      !! ** Purpose :   read wave parameters from wave model  in netcdf files.
346      !!
347      !! ** Method  : - Read namelist namsbc_wave
348      !!              - Read Cd_n10 fields in netcdf files
349      !!              - Read stokes drift 2d in netcdf files
350      !!              - Read wave number in netcdf files
351      !!              - Compute 3d stokes drift using Breivik et al.,2014
352      !!                formulation
353      !! ** action 
354      !!---------------------------------------------------------------------
355      INTEGER ::   ierror, ios   ! local integer
356      INTEGER ::   ifpr
357      !!
[9115]358      CHARACTER(len=100)     ::  cn_dir                            ! Root directory for location of drag coefficient files
[9023]359      TYPE(FLD_N), ALLOCATABLE, DIMENSION(:) ::   slf_i, slf_j     ! array of namelist informations on the fields to read
[7646]360      TYPE(FLD_N)            ::  sn_cdg, sn_usd, sn_vsd,  &
[9023]361                             &   sn_hsw, sn_wmp, sn_wfr, sn_wnum, &
[9115]362                             &   sn_tauwoc, sn_tauwx, sn_tauwy     ! informations about the fields to be read
[7646]363      !
[9023]364      NAMELIST/namsbc_wave/  sn_cdg, cn_dir, sn_usd, sn_vsd, sn_hsw, sn_wmp, sn_wfr, &
[9033]365                             sn_wnum, sn_tauwoc, sn_tauwx, sn_tauwy
[7646]366      !!---------------------------------------------------------------------
367      !
368      READ  ( numnam_ref, namsbc_wave, IOSTAT = ios, ERR = 901)
[11536]369901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namsbc_wave in reference namelist' )
[7646]370         
371      READ  ( numnam_cfg, namsbc_wave, IOSTAT = ios, ERR = 902 )
[11536]372902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namsbc_wave in configuration namelist' )
[7646]373      IF(lwm) WRITE ( numond, namsbc_wave )
374      !
375      IF( ln_cdgw ) THEN
376         IF( .NOT. cpl_wdrag ) THEN
[9115]377            ALLOCATE( sf_cd(1), STAT=ierror )               !* allocate and fill sf_wave with sn_cdg
[7646]378            IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_wave_init: unable to allocate sf_wave structure' )
[3680]379            !
380                                   ALLOCATE( sf_cd(1)%fnow(jpi,jpj,1)   )
381            IF( sn_cdg%ln_tint )   ALLOCATE( sf_cd(1)%fdta(jpi,jpj,1,2) )
[7646]382            CALL fld_fill( sf_cd, (/ sn_cdg /), cn_dir, 'sbc_wave_init', 'Wave module ', 'namsbc_wave' )
[5836]383         ENDIF
[7646]384         ALLOCATE( cdn_wave(jpi,jpj) )
385      ENDIF
386
[9033]387      IF( ln_tauwoc ) THEN
388         IF( .NOT. cpl_tauwoc ) THEN
389            ALLOCATE( sf_tauwoc(1), STAT=ierror )           !* allocate and fill sf_wave with sn_tauwoc
[7646]390            IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_wave_init: unable to allocate sf_wave structure' )
[3680]391            !
[9115]392                                     ALLOCATE( sf_tauwoc(1)%fnow(jpi,jpj,1)   )
[9033]393            IF( sn_tauwoc%ln_tint )  ALLOCATE( sf_tauwoc(1)%fdta(jpi,jpj,1,2) )
394            CALL fld_fill( sf_tauwoc, (/ sn_tauwoc /), cn_dir, 'sbc_wave_init', 'Wave module', 'namsbc_wave' )
[7646]395         ENDIF
396         ALLOCATE( tauoc_wave(jpi,jpj) )
397      ENDIF
398
[9023]399      IF( ln_tauw ) THEN
400         IF( .NOT. cpl_tauw ) THEN
401            ALLOCATE( sf_tauw(2), STAT=ierror )           !* allocate and fill sf_wave with sn_tauwx/y
402            IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_wave_init: unable to allocate sf_tauw structure' )
403            !
404            ALLOCATE( slf_j(2) )
405            slf_j(1) = sn_tauwx
406            slf_j(2) = sn_tauwy
407                                    ALLOCATE( sf_tauw(1)%fnow(jpi,jpj,1)   )
408                                    ALLOCATE( sf_tauw(2)%fnow(jpi,jpj,1)   )
409            IF( slf_j(1)%ln_tint )  ALLOCATE( sf_tauw(1)%fdta(jpi,jpj,1,2) )
410            IF( slf_j(2)%ln_tint )  ALLOCATE( sf_tauw(2)%fdta(jpi,jpj,1,2) )
411            CALL fld_fill( sf_tauw, (/ slf_j /), cn_dir, 'sbc_wave_init', 'read wave input', 'namsbc_wave' )
412         ENDIF
413         ALLOCATE( tauw_x(jpi,jpj) )
414         ALLOCATE( tauw_y(jpi,jpj) )
415      ENDIF
416
[7646]417      IF( ln_sdw ) THEN   ! Find out how many fields have to be read from file if not coupled
418         jpfld=0
[9023]419         jp_usd=0   ;   jp_vsd=0   ;   jp_hsw=0   ;   jp_wmp=0   ;   jp_wfr=0
[7646]420         IF( .NOT. cpl_sdrftx ) THEN
421            jpfld  = jpfld + 1
422            jp_usd = jpfld
423         ENDIF
424         IF( .NOT. cpl_sdrfty ) THEN
425            jpfld  = jpfld + 1
426            jp_vsd = jpfld
427         ENDIF
[9115]428         IF( .NOT. cpl_hsig  .AND. ll_st_bv_li  ) THEN
[7646]429            jpfld  = jpfld + 1
430            jp_hsw = jpfld
431         ENDIF
[9115]432         IF( .NOT. cpl_wper  .AND. ll_st_bv_li  ) THEN
[7646]433            jpfld  = jpfld + 1
434            jp_wmp = jpfld
435         ENDIF
[9115]436         IF( .NOT. cpl_wfreq .AND. ll_st_peakfr ) THEN
[9023]437            jpfld  = jpfld + 1
438            jp_wfr = jpfld
439         ENDIF
[7646]440
441         ! Read from file only the non-coupled fields
442         IF( jpfld > 0 ) THEN
443            ALLOCATE( slf_i(jpfld) )
444            IF( jp_usd > 0 )   slf_i(jp_usd) = sn_usd
445            IF( jp_vsd > 0 )   slf_i(jp_vsd) = sn_vsd
446            IF( jp_hsw > 0 )   slf_i(jp_hsw) = sn_hsw
447            IF( jp_wmp > 0 )   slf_i(jp_wmp) = sn_wmp
[9023]448            IF( jp_wfr > 0 )   slf_i(jp_wfr) = sn_wfr
449
[7646]450            ALLOCATE( sf_sd(jpfld), STAT=ierror )   !* allocate and fill sf_sd with stokes drift
451            IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_wave_init: unable to allocate sf_wave structure' )
452            !
[3680]453            DO ifpr= 1, jpfld
454               ALLOCATE( sf_sd(ifpr)%fnow(jpi,jpj,1) )
455               IF( slf_i(ifpr)%ln_tint )   ALLOCATE( sf_sd(ifpr)%fdta(jpi,jpj,1,2) )
456            END DO
[7646]457            !
458            CALL fld_fill( sf_sd, slf_i, cn_dir, 'sbc_wave_init', 'Wave module ', 'namsbc_wave' )
[3680]459         ENDIF
[7646]460         ALLOCATE( usd  (jpi,jpj,jpk), vsd  (jpi,jpj,jpk), wsd(jpi,jpj,jpk) )
461         ALLOCATE( hsw  (jpi,jpj)    , wmp  (jpi,jpj)     )
[9023]462         ALLOCATE( wfreq(jpi,jpj) )
[7646]463         ALLOCATE( ut0sd(jpi,jpj)    , vt0sd(jpi,jpj)     )
464         ALLOCATE( div_sd(jpi,jpj) )
465         ALLOCATE( tsd2d (jpi,jpj) )
[9019]466
467         ut0sd(:,:) = 0._wp
468         vt0sd(:,:) = 0._wp
469         hsw(:,:) = 0._wp
470         wmp(:,:) = 0._wp
471
[7646]472         usd(:,:,:) = 0._wp
473         vsd(:,:,:) = 0._wp
474         wsd(:,:,:) = 0._wp
[9019]475         ! Wave number needed only if ln_zdfswm=T
[7646]476         IF( .NOT. cpl_wnum ) THEN
477            ALLOCATE( sf_wn(1), STAT=ierror )           !* allocate and fill sf_wave with sn_wnum
478            IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_wave_init: unable toallocate sf_wave structure' )
479                                   ALLOCATE( sf_wn(1)%fnow(jpi,jpj,1)   )
480            IF( sn_wnum%ln_tint )  ALLOCATE( sf_wn(1)%fdta(jpi,jpj,1,2) )
481            CALL fld_fill( sf_wn, (/ sn_wnum /), cn_dir, 'sbc_wave', 'Wave module', 'namsbc_wave' )
[5836]482         ENDIF
[7646]483         ALLOCATE( wnum(jpi,jpj) )
[3680]484      ENDIF
[5836]485      !
[7646]486   END SUBROUTINE sbc_wave_init
487
[2990]488   !!======================================================================
489END MODULE sbcwave
Note: See TracBrowser for help on using the repository browser.