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 branches/UKMO/r6232_HZG_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/UKMO/r6232_HZG_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/SBC/sbcwave.F90 @ 7853

Last change on this file since 7853 was 7853, checked in by jcastill, 7 years ago

Add the original Craig and Banner vertical mixing scheme in case of wave coupling - further checks on wave forcing fields in case the ocean and wave land/sea masks are not the same

File size: 20.0 KB
Line 
1MODULE sbcwave
2   !!======================================================================
3   !!                       ***  MODULE  sbcwave  ***
4   !! Wave module
5   !!======================================================================
6   !! History :  3.3  !   2011-09  (Adani M)  Original code: Drag Coefficient
7   !!         :  3.4  !   2012-10  (Adani M)                 Stokes Drift
8   !!            3.6  !   2014-09  (Clementi E, Oddo P)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 oce            ! ocean variables
19   USE sbc_oce        ! Surface boundary condition: ocean fields
20   USE zdf_oce,  ONLY : ln_zdfqiao
21   USE bdy_oce        ! open boundary condition variables
22   USE domvvl         ! domain: variable volume layers
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   USE wrk_nemo       !
29   USE phycst         ! physical constants
30
31   IMPLICIT NONE
32   PRIVATE
33
34   PUBLIC   sbc_stokes      ! routine called in sbccpl
35   PUBLIC   sbc_wave        ! routine called in sbcmod
36   PUBLIC   sbc_wave_init   ! routine called in sbcmod
37   
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.
44   LOGICAL, PUBLIC ::   cpl_wnum   = .FALSE.
45   LOGICAL, PUBLIC ::   cpl_tauoc  = .FALSE.
46   LOGICAL, PUBLIC ::   cpl_wdrag  = .FALSE.
47
48   INTEGER ::   jpfld    ! number of files to read for stokes drift
49   INTEGER ::   jp_usd   ! index of stokes drift  (i-component) (m/s)    at T-point
50   INTEGER ::   jp_vsd   ! index of stokes drift  (j-component) (m/s)    at T-point
51   INTEGER ::   jp_hsw   ! index of significant wave hight      (m)      at T-point
52   INTEGER ::   jp_wmp   ! index of mean wave period            (s)      at T-point
53
54   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::  sf_cd    ! structure of input fields (file informations, fields read) Drag Coefficient
55   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::  sf_sd    ! structure of input fields (file informations, fields read) Stokes Drift
56   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::  sf_wn    ! structure of input fields (file informations, fields read) wave number for Qiao
57   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::  sf_tauoc ! structure of input fields (file informations, fields read) normalized wave stress into the ocean
58   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::  sf_phioc ! structure of input fields (file informations, fields read) wave to ocean energy
59   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)   ::   cdn_wave            !:
60   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)   ::   hsw, wmp, wnum      !:
61   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)   ::   rn_crban            !: Craig and Banner constant for surface breaking waves mixing
62   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)   ::   tauoc_wave          !:
63   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)   ::   tsd2d               !:
64   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)   ::   div_sd              !: barotropic stokes drift divergence
65   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)   ::   ut0sd, vt0sd        !: surface Stokes drift velocities at t-point
66   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:,:) ::   usd  , vsd  , wsd   !: Stokes drift velocities at u-, v- & w-points, resp.
67
68#  include "vectopt_loop_substitute.h90"
69   !!----------------------------------------------------------------------
70   !! NEMO/OPA 4.0 , NEMO Consortium (2011)
71   !! $Id$
72   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
73   !!----------------------------------------------------------------------
74CONTAINS
75
76   SUBROUTINE sbc_stokes( )
77      !!---------------------------------------------------------------------
78      !!                     ***  ROUTINE sbc_stokes  ***
79      !!
80      !! ** Purpose :   compute the 3d Stokes Drift according to Breivik et al.,
81      !!                2014 (DOI: 10.1175/JPO-D-14-0020.1)
82      !!
83      !! ** Method  : - Calculate Stokes transport speed
84      !!              - Calculate horizontal divergence
85      !!              - Integrate the horizontal divergenze from the bottom
86      !! ** action 
87      !!---------------------------------------------------------------------
88      INTEGER  ::   jj, ji, jk   ! dummy loop argument
89      INTEGER  ::   ik           ! local integer
90      REAL(wp) ::  ztransp, zfac, ztemp
91      REAL(wp) ::  zdep_u, zdep_v, zkh_u, zkh_v, zda_u, zda_v
92      REAL(wp), DIMENSION(:,:)  , POINTER ::   zk_t, zk_u, zk_v, zu0_sd, zv0_sd   ! 2D workspace
93      REAL(wp), DIMENSION(:,:,:), POINTER ::   ze3divh                            ! 3D workspace
94
95      !!---------------------------------------------------------------------
96      !
97
98      CALL wrk_alloc( jpi,jpj,jpk,   ze3divh )
99      CALL wrk_alloc( jpi,jpj,       zk_t, zk_u, zk_v, zu0_sd, zv0_sd ) 
100      !
101      zfac = 2.0_wp * rpi / 16.0_wp
102      DO jj = 1, jpj               ! exp. wave number at t-point    (Eq. (19) in Breivick et al. (2014) )
103         DO ji = 1, jpi
104            ! Stokes drift velocity estimated from Hs and Tmean
105            ztransp = zfac * hsw(ji,jj)*hsw(ji,jj) / MAX( wmp(ji,jj), 0.0000001_wp )
106            ! Stokes surface speed
107            tsd2d(ji,jj) = SQRT( ut0sd(ji,jj)*ut0sd(ji,jj) + vt0sd(ji,jj)*vt0sd(ji,jj))
108            ! Wavenumber scale
109            zk_t(ji,jj) = ABS( tsd2d(ji,jj) ) / MAX( ABS( 5.97_wp*ztransp ), 0.0000001_wp )
110         END DO
111      END DO
112      DO jj = 1, jpjm1              ! exp. wave number & Stokes drift velocity at u- & v-points
113         DO ji = 1, jpim1
114            zk_u(ji,jj) = 0.5_wp * ( zk_t(ji,jj) + zk_t(ji+1,jj) )
115            zk_v(ji,jj) = 0.5_wp * ( zk_t(ji,jj) + zk_t(ji,jj+1) )
116            !
117            zu0_sd(ji,jj) = 0.5_wp * ( ut0sd(ji,jj) + ut0sd(ji+1,jj) )
118            zv0_sd(ji,jj) = 0.5_wp * ( vt0sd(ji,jj) + vt0sd(ji,jj+1) )
119         END DO
120      END DO
121      !
122      !                       !==  horizontal Stokes Drift 3D velocity  ==!
123      DO jk = 1, jpkm1
124         DO jj = 2, jpjm1
125            DO ji = 2, jpim1
126               zdep_u = 0.5_wp * ( gdept_n(ji,jj,jk) + gdept_n(ji+1,jj,jk) )
127               zdep_v = 0.5_wp * ( gdept_n(ji,jj,jk) + gdept_n(ji,jj+1,jk) )
128               !                         
129               zkh_u = zk_u(ji,jj) * zdep_u     ! k * depth
130               zkh_v = zk_v(ji,jj) * zdep_v
131               !                                ! Depth attenuation
132               zda_u = EXP( -2.0_wp*zkh_u ) / ( 1.0_wp + 8.0_wp*zkh_u )
133               zda_v = EXP( -2.0_wp*zkh_v ) / ( 1.0_wp + 8.0_wp*zkh_v )
134               !
135               usd(ji,jj,jk) = zda_u * zk_u(ji,jj) * umask(ji,jj,jk)
136               vsd(ji,jj,jk) = zda_v * zk_v(ji,jj) * vmask(ji,jj,jk)
137            END DO
138         END DO
139      END DO
140      CALL lbc_lnk( usd(:,:,:), 'U', vsd(:,:,:), 'V', -1. )
141      !
142      !                       !==  vertical Stokes Drift 3D velocity  ==!
143      !
144      DO jk = 1, jpkm1               ! Horizontal e3*divergence
145         DO jj = 2, jpj
146            DO ji = fs_2, jpi
147               ze3divh(ji,jj,jk) = (  e2u(ji  ,jj) * e3u_n(ji  ,jj,jk) * usd(ji,  jj,jk)    &
148                  &                 - e2u(ji-1,jj) * e3u_n(ji-1,jj,jk) * usd(ji-1,jj,jk)    &
149                  &                 + e1v(ji,jj  ) * e3v_n(ji,jj  ,jk) * vsd(ji,jj  ,jk)    &
150                  &                 - e1v(ji,jj-1) * e3v_n(ji,jj-1,jk) * vsd(ji,jj-1,jk)  ) * r1_e12t(ji,jj)
151            END DO
152         END DO
153      END DO
154      !
155      IF( .NOT. AGRIF_Root() ) THEN
156         IF( nbondi ==  1 .OR. nbondi == 2 )   ze3divh(nlci-1,   :  ,:) = 0._wp      ! east
157         IF( nbondi == -1 .OR. nbondi == 2 )   ze3divh(  2   ,   :  ,:) = 0._wp      ! west
158         IF( nbondj ==  1 .OR. nbondj == 2 )   ze3divh(  :   ,nlcj-1,:) = 0._wp      ! north
159         IF( nbondj == -1 .OR. nbondj == 2 )   ze3divh(  :   ,  2   ,:) = 0._wp      ! south
160      ENDIF
161      !
162      CALL lbc_lnk( ze3divh, 'T', 1. )
163      !
164      IF( .NOT. lk_vvl ) THEN   ;   ik = 1   ! none zero velocity through the sea surface
165      ELSE                      ;   ik = 2   ! w=0 at the surface (set one for all in sbc_wave_init)
166      ENDIF
167      DO jk = jpkm1, ik, -1          ! integrate from the bottom the hor. divergence (NB: at k=jpk w is always zero)
168         wsd(:,:,jk) = wsd(:,:,jk+1) - ze3divh(:,:,jk)
169      END DO
170#if defined key_bdy
171      IF( lk_bdy ) THEN
172         DO jk = 1, jpkm1
173            wsd(:,:,jk) = wsd(:,:,jk) * bdytmask(:,:)
174         END DO
175      ENDIF
176#endif
177      !                       !==  Horizontal divergence of barotropic Stokes transport  ==!
178      div_sd(:,:) = 0._wp
179      DO jk = 1, jpkm1                                 !
180        div_sd(:,:) = div_sd(:,:) + ze3divh(:,:,jk)
181      END DO
182      !
183      CALL iom_put( "ustokes",  usd  )
184      CALL iom_put( "vstokes",  vsd  )
185      CALL iom_put( "wstokes",  wsd  )
186      !
187      CALL wrk_dealloc( jpi,jpj,jpk,   ze3divh )
188      CALL wrk_dealloc( jpi,jpj,       zk_t, zk_u, zk_v, zu0_sd, zv0_sd )
189      !
190   END SUBROUTINE sbc_stokes
191
192
193   SUBROUTINE sbc_wave( kt )
194      !!---------------------------------------------------------------------
195      !!                     ***  ROUTINE sbc_wave  ***
196      !!
197      !! ** Purpose :   read wave parameters from wave model  in netcdf files.
198      !!
199      !! ** Method  : - Read namelist namsbc_wave
200      !!              - Read Cd_n10 fields in netcdf files
201      !!              - Read stokes drift 2d in netcdf files
202      !!              - Read wave number in netcdf files
203      !!              - Compute 3d stokes drift using Breivik et al.,2014
204      !!                formulation
205      !! ** action 
206      !!---------------------------------------------------------------------
207      INTEGER, INTENT(in   ) ::   kt   ! ocean time step
208      !!---------------------------------------------------------------------
209      !
210      IF( ln_cdgw .AND. .NOT. cpl_wdrag ) THEN     !==  Neutral drag coefficient  ==!
211         CALL fld_read( kt, nn_fsbc, sf_cd )             ! read from external forcing
212         cdn_wave(:,:) = sf_cd(1)%fnow(:,:,1)
213         ! check that the drag coefficient contains proper information even if
214         ! the masks do not match - the momentum stress is not masked!
215         WHERE( cdn_wave < 0.0 ) cdn_wave = 1.5e-3
216         WHERE( cdn_wave > 1.0 ) cdn_wave = 1.5e-3
217      ENDIF
218
219      IF( ln_tauoc .AND. .NOT. cpl_tauoc ) THEN    !==  Wave induced stress  ==!
220         CALL fld_read( kt, nn_fsbc, sf_tauoc )          ! read wave norm stress from external forcing
221         tauoc_wave(:,:) = sf_tauoc(1)%fnow(:,:,1)
222         WHERE( tauoc_wave < -100.0 ) tauoc_wave = 1.0
223         WHERE( tauoc_wave >  100.0 ) tauoc_wave = 1.0
224      ENDIF
225
226      IF( ln_phioc .AND. .NOT. cpl_phioc ) THEN    !==  Wave to ocean energy  ==!
227         CALL fld_read( kt, nn_fsbc, sf_phioc )          ! read wave to ocean energy from external forcing
228         rn_crban(:,:) = 29.0 * sf_phioc(1)%fnow(:,:,1)     ! ! Alfa is phioc*sqrt(rau0/zrhoa)  : rau0=water density, zhroa= air density
229         WHERE( rn_crban <  10.0 ) rn_crban =  10.0
230         WHERE( rn_crban > 300.0 ) rn_crban = 300.0
231      ENDIF
232
233      IF( ln_sdw )  THEN                           !==  Computation of the 3d Stokes Drift  ==!
234         !
235         IF( jpfld > 0 ) THEN                            ! Read from file only if the field is not coupled
236            CALL fld_read( kt, nn_fsbc, sf_sd )          ! read wave parameters from external forcing
237            IF( jp_hsw > 0 ) THEN
238               hsw  (:,:) = sf_sd(jp_hsw)%fnow(:,:,1)   ! significant wave height
239               WHERE( hsw > 100.0 ) hsw = 0.0
240               WHERE( hsw <   0.0 ) hsw = 0.0
241            ENDIF
242            IF( jp_wmp > 0 ) THEN
243               wmp  (:,:) = sf_sd(jp_wmp)%fnow(:,:,1)   ! wave mean period
244               WHERE( wmp > 100.0 ) wmp = 0.0
245               WHERE( wmp <   0.0 ) wmp = 0.0
246            ENDIF
247            IF( jp_usd > 0 ) THEN
248               ut0sd(:,:) = sf_sd(jp_usd)%fnow(:,:,1)   ! 2D zonal Stokes Drift at T point
249               WHERE( ut0sd < -100.0 ) ut0sd = 1.0
250               WHERE( ut0sd >  100.0 ) ut0sd = 1.0
251            ENDIF
252            IF( jp_vsd > 0 ) THEN
253               vt0sd(:,:) = sf_sd(jp_vsd)%fnow(:,:,1)   ! 2D meridional Stokes Drift at T point
254               WHERE( vt0sd < -100.0 ) vt0sd = 1.0
255               WHERE( vt0sd >  100.0 ) vt0sd = 1.0
256            ENDIF
257         ENDIF
258         !
259         ! Read also wave number if needed, so that it is available in coupling routines
260         IF( ln_zdfqiao .AND. .NOT.cpl_wnum ) THEN
261            CALL fld_read( kt, nn_fsbc, sf_wn )          ! read wave parameters from external forcing
262            wnum(:,:) = sf_wn(1)%fnow(:,:,1)
263         ENDIF
264           
265         !                                         !==  Computation of the 3d Stokes Drift  ==!
266         !
267         IF( jpfld == 4 )   CALL sbc_stokes()            ! Calculate only if required fields are read
268         !                                               ! In coupled wave model-NEMO case the call is done after coupling
269         !
270      ENDIF
271      !
272   END SUBROUTINE sbc_wave
273
274
275   SUBROUTINE sbc_wave_init
276      !!---------------------------------------------------------------------
277      !!                     ***  ROUTINE sbc_wave_init  ***
278      !!
279      !! ** Purpose :   read wave parameters from wave model  in netcdf files.
280      !!
281      !! ** Method  : - Read namelist namsbc_wave
282      !!              - Read Cd_n10 fields in netcdf files
283      !!              - Read stokes drift 2d in netcdf files
284      !!              - Read wave number in netcdf files
285      !!              - Compute 3d stokes drift using Breivik et al.,2014
286      !!                formulation
287      !! ** action 
288      !!---------------------------------------------------------------------
289      INTEGER ::   ierror, ios   ! local integer
290      INTEGER ::   ifpr
291      !!
292      CHARACTER(len=100)     ::  cn_dir                          ! Root directory for location of drag coefficient files
293      TYPE(FLD_N), ALLOCATABLE, DIMENSION(:) ::   slf_i     ! array of namelist informations on the fields to read
294      TYPE(FLD_N)            ::  sn_cdg, sn_usd, sn_vsd,  sn_phioc, &
295                             &   sn_hsw, sn_wmp, sn_wnum, sn_tauoc      ! informations about the fields to be read
296      !
297      NAMELIST/namsbc_wave/  sn_cdg, cn_dir, sn_usd, sn_vsd, sn_hsw, sn_wmp, sn_wnum, sn_tauoc, sn_phioc
298      !!---------------------------------------------------------------------
299      !
300      REWIND( numnam_ref )              ! Namelist namsbc_wave in reference namelist : File for drag coeff. from wave model
301      READ  ( numnam_ref, namsbc_wave, IOSTAT = ios, ERR = 901)
302901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_wave in reference namelist', lwp )
303         
304      REWIND( numnam_cfg )              ! Namelist namsbc_wave in configuration namelist : File for drag coeff. from wave model
305      READ  ( numnam_cfg, namsbc_wave, IOSTAT = ios, ERR = 902 )
306902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_wave in configuration namelist', lwp )
307      IF(lwm) WRITE ( numond, namsbc_wave )
308      !
309      IF( ln_cdgw ) THEN
310         IF( .NOT. cpl_wdrag ) THEN
311            ALLOCATE( sf_cd(1), STAT=ierror )           !* allocate and fill sf_wave with sn_cdg
312            IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_wave_init: unable to allocate sf_cd structure' )
313            !
314                                   ALLOCATE( sf_cd(1)%fnow(jpi,jpj,1)   )
315            IF( sn_cdg%ln_tint )   ALLOCATE( sf_cd(1)%fdta(jpi,jpj,1,2) )
316            CALL fld_fill( sf_cd, (/ sn_cdg /), cn_dir, 'sbc_wave_init', 'Wave module ', 'namsbc_wave' )
317         ENDIF
318         ALLOCATE( cdn_wave(jpi,jpj) )
319      ENDIF
320
321      IF( ln_tauoc ) THEN
322         IF( .NOT. cpl_tauoc ) THEN
323            ALLOCATE( sf_tauoc(1), STAT=ierror )           !* allocate and fill sf_wave with sn_tauoc
324            IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_wave_init: unable to allocate sf_tauoc structure' )
325            !
326                                    ALLOCATE( sf_tauoc(1)%fnow(jpi,jpj,1)   )
327            IF( sn_tauoc%ln_tint )  ALLOCATE( sf_tauoc(1)%fdta(jpi,jpj,1,2) )
328            CALL fld_fill( sf_tauoc, (/ sn_tauoc /), cn_dir, 'sbc_wave_init', 'Wave module', 'namsbc_wave' )
329         ENDIF
330         ALLOCATE( tauoc_wave(jpi,jpj) )
331      ENDIF
332
333      IF( ln_phioc ) THEN
334         IF( .NOT. cpl_phioc ) THEN
335            ALLOCATE( sf_phioc(1), STAT=ierror )           !* allocate and fill sf_wave with sn_phioc
336            IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_wave_init: unable to allocate sf_phioc structure' )
337            !
338                                    ALLOCATE( sf_phioc(1)%fnow(jpi,jpj,1)   )
339            IF( sn_phioc%ln_tint )  ALLOCATE( sf_phioc(1)%fdta(jpi,jpj,1,2) )
340            CALL fld_fill( sf_phioc, (/ sn_phioc /), cn_dir, 'sbc_wave_init', 'Wave module', 'namsbc_wave' )
341         ENDIF
342         ALLOCATE( rn_crban(jpi,jpj) )
343      ENDIF
344
345      IF( ln_sdw ) THEN   ! Find out how many fields have to be read from file if not coupled
346         jpfld=0
347         jp_usd=0   ;   jp_vsd=0   ;   jp_hsw=0   ;   jp_wmp=0
348         IF( .NOT. cpl_sdrftx ) THEN
349            jpfld  = jpfld + 1
350            jp_usd = jpfld
351         ENDIF
352         IF( .NOT. cpl_sdrfty ) THEN
353            jpfld  = jpfld + 1
354            jp_vsd = jpfld
355         ENDIF
356         IF( .NOT. cpl_hsig ) THEN
357            jpfld  = jpfld + 1
358            jp_hsw = jpfld
359         ENDIF
360         IF( .NOT. cpl_wper ) THEN
361            jpfld  = jpfld + 1
362            jp_wmp = jpfld
363         ENDIF
364
365         ! Read from file only the non-coupled fields
366         IF( jpfld > 0 ) THEN
367            ALLOCATE( slf_i(jpfld) )
368            IF( jp_usd > 0 )   slf_i(jp_usd) = sn_usd
369            IF( jp_vsd > 0 )   slf_i(jp_vsd) = sn_vsd
370            IF( jp_hsw > 0 )   slf_i(jp_hsw) = sn_hsw
371            IF( jp_wmp > 0 )   slf_i(jp_wmp) = sn_wmp
372            ALLOCATE( sf_sd(jpfld), STAT=ierror )   !* allocate and fill sf_sd with stokes drift
373            IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_wave_init: unable to allocate sf_sd structure' )
374            !
375            DO ifpr= 1, jpfld
376               ALLOCATE( sf_sd(ifpr)%fnow(jpi,jpj,1) )
377               IF( slf_i(ifpr)%ln_tint )   ALLOCATE( sf_sd(ifpr)%fdta(jpi,jpj,1,2) )
378            END DO
379            !
380            CALL fld_fill( sf_sd, slf_i, cn_dir, 'sbc_wave_init', 'Wave module ', 'namsbc_wave' )
381         ENDIF
382         ALLOCATE( usd  (jpi,jpj,jpk), vsd  (jpi,jpj,jpk), wsd(jpi,jpj,jpk) )
383         ALLOCATE( hsw  (jpi,jpj)    , wmp  (jpi,jpj)     )
384         ALLOCATE( ut0sd(jpi,jpj)    , vt0sd(jpi,jpj)     )
385         ALLOCATE( div_sd(jpi,jpj) )
386         ALLOCATE( tsd2d (jpi,jpj) )
387         usd(:,:,:) = 0._wp
388         vsd(:,:,:) = 0._wp
389         wsd(:,:,:) = 0._wp
390         ! Wave number needed only if ln_zdfqiao=T
391         IF( ln_zdfqiao .AND. .NOT.cpl_wnum ) THEN
392            ALLOCATE( sf_wn(1), STAT=ierror )           !* allocate and fill sf_wave with sn_wnum
393            IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_wave_init: unable toallocate sf_wn structure' )
394                                   ALLOCATE( sf_wn(1)%fnow(jpi,jpj,1)   )
395            IF( sn_wnum%ln_tint )  ALLOCATE( sf_wn(1)%fdta(jpi,jpj,1,2) )
396            CALL fld_fill( sf_wn, (/ sn_wnum /), cn_dir, 'sbc_wave', 'Wave module', 'namsbc_wave' )
397         ENDIF
398         ALLOCATE( wnum(jpi,jpj) )
399      ENDIF
400      !
401   END SUBROUTINE sbc_wave_init
402
403   !!======================================================================
404END MODULE sbcwave
Note: See TracBrowser for help on using the repository browser.