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/2018/dev_r10057_ENHANCE03_ZTILDE/src/OCE/SBC – NEMO

source: NEMO/branches/2018/dev_r10057_ENHANCE03_ZTILDE/src/OCE/SBC/sbcwave.F90 @ 10116

Last change on this file since 10116 was 10116, checked in by jchanut, 6 years ago

ztilde update (3): #2126
Rewrite thickness advection, add regriding, biharmonic interface diffusion

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