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.
zdfbfr.F90 in branches/2014/dev_r4650_UKMO2_ice_shelves/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

source: branches/2014/dev_r4650_UKMO2_ice_shelves/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfbfr.F90 @ 4724

Last change on this file since 4724 was 4724, checked in by mathiot, 10 years ago

ISF branch: add comments, fix mpp and restar issues, add test to stop if incompatible options and fix mask issue in sbcice and sbcblk.

  • Property svn:keywords set to Id
File size: 23.1 KB
Line 
1MODULE zdfbfr
2   !!======================================================================
3   !!                       ***  MODULE  zdfbfr  ***
4   !! Ocean physics: Bottom friction
5   !!======================================================================
6   !! History :  OPA  ! 1997-06  (G. Madec, A.-M. Treguier)  Original code
7   !!   NEMO     1.0  ! 2002-06  (G. Madec)  F90: Free form and module
8   !!            3.2  ! 2009-09  (A.C.Coward)  Correction to include barotropic contribution
9   !!            3.3  ! 2010-10  (C. Ethe, G. Madec) reorganisation of initialisation phase
10   !!            3.4  ! 2011-11  (H. Liu) implementation of semi-implicit bottom friction option
11   !!                 ! 2012-06  (H. Liu) implementation of Log Layer bottom friction option
12   !!----------------------------------------------------------------------
13
14   !!----------------------------------------------------------------------
15   !!   zdf_bfr      : update bottom friction coefficient (non-linear bottom friction only)
16   !!   zdf_bfr_init : read in namelist and control the bottom friction parameters.
17   !!----------------------------------------------------------------------
18   USE oce             ! ocean dynamics and tracers variables
19   USE dom_oce         ! ocean space and time domain variables
20   USE zdf_oce         ! ocean vertical physics variables
21   USE in_out_manager  ! I/O manager
22   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
23   USE lib_mpp         ! distributed memory computing
24   USE prtctl          ! Print control
25   USE timing          ! Timing
26   USE wrk_nemo        ! Memory Allocation
27   USE phycst, ONLY: vkarmn
28
29   IMPLICIT NONE
30   PRIVATE
31
32   PUBLIC   zdf_bfr         ! called by step.F90
33   PUBLIC   zdf_bfr_init    ! called by nemogcm.F90
34
35   !                                 !!* Namelist nambfr: bottom friction namelist *
36   INTEGER , PUBLIC ::   nn_bfr       ! = 0/1/2/3 type of bottom friction  (PUBLIC for TAM)
37   REAL(wp), PUBLIC ::   rn_bfri1     ! bottom drag coefficient (linear case)  (PUBLIC for TAM)
38   REAL(wp), PUBLIC ::   rn_bfri2     ! bottom drag coefficient (non linear case) (PUBLIC for TAM)
39   REAL(wp), PUBLIC ::   rn_bfri2_max ! Maximum bottom drag coefficient (non linear case and ln_loglayer=T) (PUBLIC for TAM)
40   REAL(wp), PUBLIC ::   rn_bfeb2     ! background bottom turbulent kinetic energy  [m2/s2] (PUBLIC for TAM)
41   REAL(wp), PUBLIC ::   rn_bfrien    ! local factor to enhance coefficient bfri (PUBLIC for TAM)
42   LOGICAL , PUBLIC ::   ln_bfr2d     ! logical switch for 2D enhancement (PUBLIC for TAM)
43   REAL(wp), PUBLIC ::   rn_tfri1     ! top drag coefficient (linear case)  (PUBLIC for TAM)
44   REAL(wp), PUBLIC ::   rn_tfri2     ! top drag coefficient (non linear case) (PUBLIC for TAM)
45   REAL(wp), PUBLIC ::   rn_tfri2_max ! Maximum top drag coefficient (non linear case and ln_loglayer=T) (PUBLIC for TAM)
46   REAL(wp), PUBLIC ::   rn_tfeb2     ! background top turbulent kinetic energy  [m2/s2] (PUBLIC for TAM)
47   REAL(wp), PUBLIC ::   rn_tfrien    ! local factor to enhance coefficient tfri (PUBLIC for TAM)
48   LOGICAL , PUBLIC ::   ln_tfr2d     ! logical switch for 2D enhancement (PUBLIC for TAM)
49
50   LOGICAL , PUBLIC ::   ln_loglayer  ! switch for log layer bfr coeff. (PUBLIC for TAM)
51   REAL(wp), PUBLIC ::   rn_bfrz0     ! bottom roughness for loglayer bfr coeff (PUBLIC for TAM)
52   REAL(wp), PUBLIC ::   rn_tfrz0     ! bottom roughness for loglayer bfr coeff (PUBLIC for TAM)
53   LOGICAL , PUBLIC ::   ln_bfrimp    ! logical switch for implicit bottom friction
54   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:), PUBLIC ::  bfrcoef2d, tfrcoef2d   ! 2D bottom/top drag coefficient (PUBLIC for TAM)
55
56   !! * Substitutions
57#  include "vectopt_loop_substitute.h90"
58#  include "domzgr_substitute.h90"
59   !!----------------------------------------------------------------------
60   !! NEMO/OPA 4.0 , NEMO Consortium (2011)
61   !! $Id$
62   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
63   !!----------------------------------------------------------------------
64CONTAINS
65
66   INTEGER FUNCTION zdf_bfr_alloc()
67      !!----------------------------------------------------------------------
68      !!                ***  FUNCTION zdf_bfr_alloc  ***
69      !!----------------------------------------------------------------------
70      ALLOCATE( bfrcoef2d(jpi,jpj), tfrcoef2d(jpi,jpj), STAT=zdf_bfr_alloc )
71      !
72      IF( lk_mpp             )   CALL mpp_sum ( zdf_bfr_alloc )
73      IF( zdf_bfr_alloc /= 0 )   CALL ctl_warn('zdf_bfr_alloc: failed to allocate arrays.')
74   END FUNCTION zdf_bfr_alloc
75
76
77   SUBROUTINE zdf_bfr( kt )
78      !!----------------------------------------------------------------------
79      !!                   ***  ROUTINE zdf_bfr  ***
80      !!
81      !! ** Purpose :   compute the bottom friction coefficient.
82      !!
83      !! ** Method  :   Calculate and store part of the momentum trend due
84      !!              to bottom friction following the chosen friction type
85      !!              (free-slip, linear, or quadratic). The component
86      !!              calculated here is multiplied by the bottom velocity in
87      !!              dyn_bfr to provide the trend term.
88      !!                The coefficients are updated at each time step only
89      !!              in the quadratic case.
90      !!
91      !! ** Action  :   bfrua , bfrva   bottom friction coefficients
92      !!----------------------------------------------------------------------
93      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index
94      !!
95      INTEGER  ::   ji, jj                       ! dummy loop indices
96      INTEGER  ::   ikbt, ikbu, ikbv             ! local integers
97      REAL(wp) ::   zvu, zuv, zecu, zecv, ztmp   ! temporary scalars
98      REAL(wp), POINTER, DIMENSION(:,:) ::  zbfrt, ztfrt
99      !!----------------------------------------------------------------------
100      !
101      IF( nn_timing == 1 )  CALL timing_start('zdf_bfr')
102      !
103      IF( kt == nit000 .AND. lwp ) THEN
104         WRITE(numout,*)
105         WRITE(numout,*) 'zdf_bfr : Set bottom friction coefficient (non-linear case)'
106         WRITE(numout,*) '~~~~~~~~'
107      ENDIF
108      !
109      IF( nn_bfr == 2 ) THEN                 ! quadratic bottom friction only
110         !
111         CALL wrk_alloc( jpi, jpj, zbfrt, ztfrt )
112
113         IF ( ln_loglayer.AND.lk_vvl ) THEN ! "log layer" bottom friction coefficient
114
115#  if defined key_vectopt_loop
116            DO jj = 1, 1
117!CDIR NOVERRCHK
118               DO ji = 1, jpij   ! vector opt. (forced unrolling)
119#  else
120!CDIR NOVERRCHK
121            DO jj = 1, jpj
122!CDIR NOVERRCHK
123               DO ji = 1, jpi
124#  endif
125                  ikbt = mbkt(ji,jj)
126! JC: possible WAD implementation should modify line below if layers vanish
127                  ztmp = tmask(ji,jj,ikbt) * ( vkarmn / LOG( 0.5_wp * fse3t_n(ji,jj,ikbt) / rn_bfrz0 ))**2._wp
128                  zbfrt(ji,jj) = MAX(bfrcoef2d(ji,jj), ztmp)
129                  zbfrt(ji,jj) = MIN(zbfrt(ji,jj), rn_bfri2_max)
130! (ISF)
131                  ikbt = mikt(ji,jj)
132! JC: possible WAD implementation should modify line below if layers vanish
133                  ztmp = tmask(ji,jj,ikbt) * ( vkarmn / LOG( 0.5_wp * fse3t_n(ji,jj,ikbt) / rn_bfrz0 ))**2._wp
134                  ztfrt(ji,jj) = MAX(tfrcoef2d(ji,jj), ztmp)
135                  ztfrt(ji,jj) = MIN(ztfrt(ji,jj), rn_tfri2_max)
136
137               END DO
138            END DO
139         !   
140         ELSE
141            zbfrt(:,:) = bfrcoef2d(:,:)
142            ztfrt(:,:) = tfrcoef2d(:,:)
143         ENDIF
144
145# if defined key_vectopt_loop
146         DO jj = 1, 1
147!CDIR NOVERRCHK
148            DO ji = jpi+2, jpij-jpi-1   ! vector opt. (forced unrolling)
149# else
150!CDIR NOVERRCHK
151         DO jj = 2, jpjm1
152!CDIR NOVERRCHK
153            DO ji = 2, jpim1
154# endif
155               ikbu = mbku(ji,jj)         ! ocean bottom level at u- and v-points
156               ikbv = mbkv(ji,jj)         ! (deepest ocean u- and v-points)
157               !
158               zvu  = 0.25 * (  vn(ji,jj  ,ikbu) + vn(ji+1,jj  ,ikbu)     &
159                  &           + vn(ji,jj-1,ikbu) + vn(ji+1,jj-1,ikbu)  )
160               zuv  = 0.25 * (  un(ji,jj  ,ikbv) + un(ji-1,jj  ,ikbv)     &
161                  &           + un(ji,jj+1,ikbv) + un(ji-1,jj+1,ikbv)  )
162               !
163               zecu = SQRT(  un(ji,jj,ikbu) * un(ji,jj,ikbu) + zvu*zvu + rn_bfeb2  )
164               zecv = SQRT(  vn(ji,jj,ikbv) * vn(ji,jj,ikbv) + zuv*zuv + rn_bfeb2  )
165               !
166               bfrua(ji,jj) = - 0.5_wp * ( zbfrt(ji,jj) + zbfrt(ji+1,jj  ) ) * zecu
167               bfrva(ji,jj) = - 0.5_wp * ( zbfrt(ji,jj) + zbfrt(ji  ,jj+1) ) * zecv
168               !
169               ! (ISF) ========================================================================
170               ikbu = miku(ji,jj)         ! ocean bottom level at u- and v-points
171               ikbv = mikv(ji,jj)         ! (deepest ocean u- and v-points)
172               !
173               zvu  = 0.25 * (  vn(ji,jj  ,ikbu) + vn(ji+1,jj  ,ikbu)     &
174                  &           + vn(ji,jj-1,ikbu) + vn(ji+1,jj-1,ikbu)  )
175               zuv  = 0.25 * (  un(ji,jj  ,ikbv) + un(ji-1,jj  ,ikbv)     &
176                  &           + un(ji,jj+1,ikbv) + un(ji-1,jj+1,ikbv)  )
177               !
178               zecu = SQRT(  un(ji,jj,ikbu) * un(ji,jj,ikbu) + zvu*zvu + rn_bfeb2 )
179               zecv = SQRT(  vn(ji,jj,ikbv) * vn(ji,jj,ikbv) + zuv*zuv + rn_bfeb2 )
180               !
181               tfrua(ji,jj) = - 0.5_wp * ( ztfrt(ji,jj) + ztfrt(ji+1,jj  ) ) * zecu * (1._wp - umask(ji,jj,1))
182               tfrva(ji,jj) = - 0.5_wp * ( ztfrt(ji,jj) + ztfrt(ji  ,jj+1) ) * zecv * (1._wp - vmask(ji,jj,1))
183               ! (ISF) END ====================================================================
184
185            END DO
186         END DO
187
188         !
189         CALL lbc_lnk( bfrua, 'U', 1. )   ;   CALL lbc_lnk( bfrva, 'V', 1. )      ! Lateral boundary condition
190         CALL lbc_lnk( tfrua, 'U', 1. )   ;   CALL lbc_lnk( tfrva, 'V', 1. )      ! Lateral boundary condition
191         !
192         IF(ln_ctl)   CALL prt_ctl( tab2d_1=bfrua, clinfo1=' bfr  - u: ', mask1=umask,        &
193            &                       tab2d_2=bfrva, clinfo2=       ' v: ', mask2=vmask,ovlap=1 )
194         CALL wrk_dealloc( jpi,jpj,zbfrt, ztfrt )
195      ENDIF
196      !
197      IF( nn_timing == 1 )  CALL timing_stop('zdf_bfr')
198      !
199   END SUBROUTINE zdf_bfr
200
201
202   SUBROUTINE zdf_bfr_init
203      !!----------------------------------------------------------------------
204      !!                  ***  ROUTINE zdf_bfr_init  ***
205      !!
206      !! ** Purpose :   Initialization of the bottom friction
207      !!
208      !! ** Method  :   Read the nambfr namelist and check their consistency
209      !!                called at the first timestep (nit000)
210      !!----------------------------------------------------------------------
211      USE iom   ! I/O module for ehanced bottom friction file
212      !!
213      INTEGER   ::   inum         ! logical unit for enhanced bottom friction file
214      INTEGER   ::   ji, jj       ! dummy loop indexes
215      INTEGER   ::   ikbt, ikbu, ikbv   ! temporary integers
216      INTEGER   ::   ictu, ictv         !    -          -
217      INTEGER   ::   ios
218      REAL(wp)  ::   zminbfr, zmaxbfr   ! temporary scalars
219      REAL(wp)  ::   zmintfr, zmaxtfr   ! temporary scalars
220      REAL(wp)  ::   ztmp, zfru, zfrv   !    -         -
221      !!
222      NAMELIST/nambfr/ nn_bfr, rn_bfri1, rn_bfri2, rn_bfri2_max, rn_bfeb2, rn_bfrz0, ln_bfr2d, &
223                    &          rn_tfri1, rn_tfri2, rn_tfri2_max, rn_tfeb2, rn_tfrz0, ln_tfr2d, &
224                    &  rn_bfrien, rn_tfrien, ln_bfrimp, ln_loglayer
225      !!----------------------------------------------------------------------
226      !
227      IF( nn_timing == 1 )  CALL timing_start('zdf_bfr_init')
228      !
229      !                              !* Allocate zdfbfr arrays
230      IF( zdf_bfr_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_bfr_init : unable to allocate arrays' )
231      !
232      !                              !* Parameter control and print
233      !
234      REWIND( numnam_ref )              ! Namelist nambfr in reference namelist : Bottom momentum boundary condition
235      READ  ( numnam_ref, nambfr, IOSTAT = ios, ERR = 901)
236901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nambfr in reference namelist', lwp )
237
238      REWIND( numnam_cfg )              ! Namelist nambfr in configuration namelist : Bottom momentum boundary condition
239      READ  ( numnam_cfg, nambfr, IOSTAT = ios, ERR = 902 )
240902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nambfr in configuration namelist', lwp )
241      IF(lwm) WRITE ( numond, nambfr )
242      IF(lwp) WRITE(numout,*)
243      IF(lwp) WRITE(numout,*) 'zdf_bfr_init : momentum bottom friction'
244      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~'
245      IF(lwp) WRITE(numout,*) '   Namelist nam_bfr : set bottom friction parameters'
246      !
247      SELECT CASE (nn_bfr)
248      !
249      CASE( 0 )
250         IF(lwp) WRITE(numout,*) '      free-slip '
251         bfrua(:,:) = 0.e0
252         bfrva(:,:) = 0.e0
253         tfrua(:,:) = 0.e0
254         tfrva(:,:) = 0.e0
255         !
256      CASE( 1 )
257         IF(lwp) WRITE(numout,*) '      linear botton friction'
258         IF(lwp) WRITE(numout,*) '      bottom friction coef.   rn_bfri1  = ', rn_bfri1
259         IF( ln_bfr2d ) THEN
260            IF(lwp) WRITE(numout,*) '      read coef. enhancement distribution from file   ln_bfr2d  = ', ln_bfr2d
261            IF(lwp) WRITE(numout,*) '      coef rn_bfri2 enhancement factor                rn_bfrien  = ',rn_bfrien
262         ENDIF
263         IF(lwp) WRITE(numout,*) '      top    friction coef.   rn_bfri1  = ', rn_bfri1
264         IF( ln_tfr2d ) THEN
265            IF(lwp) WRITE(numout,*) '      read coef. enhancement distribution from file   ln_tfr2d  = ', ln_tfr2d
266            IF(lwp) WRITE(numout,*) '      coef rn_tfri2 enhancement factor                rn_tfrien  = ',rn_tfrien
267         ENDIF
268         !
269         IF(ln_bfr2d) THEN
270            ! bfr_coef is a coefficient in [0,1] giving the mask where to apply the bfr enhancement
271            CALL iom_open('bfr_coef.nc',inum)
272            CALL iom_get (inum, jpdom_data, 'bfr_coef',bfrcoef2d,1) ! bfrcoef2d is used as tmp array
273            CALL iom_close(inum)
274            bfrcoef2d(:,:) = rn_bfri1 * ( 1 + rn_bfrien * bfrcoef2d(:,:) )
275         ELSE
276            bfrcoef2d(:,:) = rn_bfri1  ! initialize bfrcoef2d to the namelist variable
277         ENDIF
278         !
279         bfrua(:,:) = - bfrcoef2d(:,:)
280         bfrva(:,:) = - bfrcoef2d(:,:)
281         !
282         IF(ln_tfr2d) THEN
283            ! tfr_coef is a coefficient in [0,1] giving the mask where to apply the bfr enhancement
284            CALL iom_open('tfr_coef.nc',inum)
285            CALL iom_get (inum, jpdom_data, 'tfr_coef',tfrcoef2d,1) ! tfrcoef2d is used as tmp array
286            CALL iom_close(inum)
287            tfrcoef2d(:,:) = rn_tfri1 * ( 1 + rn_tfrien * tfrcoef2d(:,:) )
288         ELSE
289            tfrcoef2d(:,:) = rn_tfri1  ! initialize tfrcoef2d to the namelist variable
290         ENDIF
291         !
292         tfrua(:,:) = - tfrcoef2d(:,:)
293         tfrva(:,:) = - tfrcoef2d(:,:)
294         !
295      CASE( 2 )
296         IF(lwp) WRITE(numout,*) '      quadratic bottom friction'
297         IF(lwp) WRITE(numout,*) '      friction coef.   rn_bfri2  = ', rn_bfri2
298         IF(lwp) WRITE(numout,*) '      Max. coef. (log case)   rn_bfri2_max  = ', rn_bfri2_max
299         IF(lwp) WRITE(numout,*) '      background tke   rn_bfeb2  = ', rn_bfeb2
300         IF(lwp) WRITE(numout,*) '      log formulation   ln_bfr2d = ', ln_loglayer
301         IF(lwp) WRITE(numout,*) '      bottom roughness  rn_bfrz0 [m] = ', rn_bfrz0
302         IF( rn_bfrz0<=0.e0 ) THEN
303            WRITE(ctmp1,*) '      bottom roughness must be strictly positive'
304            CALL ctl_stop( ctmp1 )
305         ENDIF
306         IF( ln_bfr2d ) THEN
307            IF(lwp) WRITE(numout,*) '      read coef. enhancement distribution from file   ln_bfr2d  = ', ln_bfr2d
308            IF(lwp) WRITE(numout,*) '      coef rn_bfri2 enhancement factor                rn_bfrien  = ',rn_bfrien
309         ENDIF
310         IF(lwp) WRITE(numout,*) '      quadratic top    friction'
311         IF(lwp) WRITE(numout,*) '      friction coef.   rn_bfri2  = ', rn_tfri2
312         IF(lwp) WRITE(numout,*) '      Max. coef. (log case)   rn_tfri2_max  = ', rn_tfri2_max
313         IF(lwp) WRITE(numout,*) '      background tke   rn_tfeb2  = ', rn_tfeb2
314         IF(lwp) WRITE(numout,*) '      log formulation   ln_tfr2d = ', ln_loglayer
315         IF(lwp) WRITE(numout,*) '      bottom roughness  rn_tfrz0 [m] = ', rn_tfrz0
316         IF( rn_tfrz0<=0.e0 ) THEN
317            WRITE(ctmp1,*) '      bottom roughness must be strictly positive'
318            CALL ctl_stop( ctmp1 )
319         ENDIF
320         IF( ln_tfr2d ) THEN
321            IF(lwp) WRITE(numout,*) '      read coef. enhancement distribution from file   ln_tfr2d  = ', ln_tfr2d
322            IF(lwp) WRITE(numout,*) '      coef rn_tfri2 enhancement factor                rn_tfrien  = ',rn_tfrien
323         ENDIF
324         !
325         IF(ln_bfr2d) THEN
326            ! bfr_coef is a coefficient in [0,1] giving the mask where to apply the bfr enhancement
327            CALL iom_open('bfr_coef.nc',inum)
328            CALL iom_get (inum, jpdom_data, 'bfr_coef',bfrcoef2d,1) ! bfrcoef2d is used as tmp array
329            CALL iom_close(inum)
330            !
331            bfrcoef2d(:,:) = rn_bfri2 * ( 1 + rn_bfrien * bfrcoef2d(:,:) )
332         ELSE
333            bfrcoef2d(:,:) = rn_bfri2  ! initialize bfrcoef2d to the namelist variable
334         ENDIF
335
336         IF(ln_tfr2d) THEN
337            ! tfr_coef is a coefficient in [0,1] giving the mask where to apply the bfr enhancement
338            CALL iom_open('tfr_coef.nc',inum)
339            CALL iom_get (inum, jpdom_data, 'tfr_coef',tfrcoef2d,1) ! bfrcoef2d is used as tmp array
340            CALL iom_close(inum)
341            !
342            tfrcoef2d(:,:) = rn_tfri2 * ( 1 + rn_tfrien * tfrcoef2d(:,:) )
343         ELSE
344            tfrcoef2d(:,:) = rn_tfri2  ! initialize tfrcoef2d to the namelist variable
345         ENDIF
346         !
347         IF ( ln_loglayer.AND.(.NOT.lk_vvl) ) THEN ! set "log layer" bottom friction once for all
348#  if defined key_vectopt_loop
349            DO jj = 1, 1
350!CDIR NOVERRCHK
351               DO ji = 1, jpij   ! vector opt. (forced unrolling)
352#  else
353!CDIR NOVERRCHK
354            DO jj = 1, jpj
355!CDIR NOVERRCHK
356               DO ji = 1, jpi
357#  endif
358                  ikbt = mbkt(ji,jj)
359                  ztmp = tmask(ji,jj,ikbt) * ( vkarmn / LOG( 0.5_wp * fse3t_n(ji,jj,ikbt) / rn_bfrz0 ))**2._wp
360                  bfrcoef2d(ji,jj) = MAX(bfrcoef2d(ji,jj), ztmp)
361                  bfrcoef2d(ji,jj) = MIN(bfrcoef2d(ji,jj), rn_bfri2_max)
362                  ikbt = mikt(ji,jj)
363                  ztmp = tmask(ji,jj,ikbt) * ( vkarmn / LOG( 0.5_wp * fse3t_n(ji,jj,ikbt) / rn_tfrz0 ))**2._wp
364                  tfrcoef2d(ji,jj) = MAX(tfrcoef2d(ji,jj), ztmp)
365                  tfrcoef2d(ji,jj) = MIN(tfrcoef2d(ji,jj), rn_tfri2_max)
366               END DO
367            END DO
368         ENDIF
369         !
370      CASE DEFAULT
371         IF(lwp) WRITE(ctmp1,*) '         bad flag value for nn_bfr = ', nn_bfr
372         CALL ctl_stop( ctmp1 )
373         !
374      END SELECT
375      !
376      IF(lwp) WRITE(numout,*) '      implicit bottom friction switch                ln_bfrimp  = ', ln_bfrimp
377      !
378      !                              ! Make sure ln_zdfexp=.false. when use implicit bfr
379      IF( ln_bfrimp .AND. ln_zdfexp ) THEN
380         IF(lwp) THEN
381            WRITE(numout,*)
382            WRITE(numout,*) 'Implicit bottom friction can only be used when ln_zdfexp=.false.'
383            WRITE(numout,*) '         but you set: ln_bfrimp=.true. and ln_zdfexp=.true.!!!!'
384            WRITE(ctmp1,*)  '         set either ln_zdfexp = .false or ln_bfrimp = .false.'
385            CALL ctl_stop( ctmp1 )
386         END IF
387      END IF
388      !
389      ! Basic stability check on bottom friction coefficient
390      !
391      ictu = 0               ! counter for stability criterion breaches at U-pts
392      ictv = 0               ! counter for stability criterion breaches at V-pts
393      zminbfr =  1.e10_wp    ! initialise tracker for minimum of bottom friction coefficient
394      zmaxbfr = -1.e10_wp    ! initialise tracker for maximum of bottom friction coefficient
395      zmintfr =  1.e10_wp    ! initialise tracker for minimum of bottom friction coefficient
396      zmaxtfr = -1.e10_wp    ! initialise tracker for maximum of bottom friction coefficient
397      !
398#  if defined key_vectopt_loop
399      DO jj = 1, 1
400!CDIR NOVERRCHK
401         DO ji = jpi+2, jpij-jpi-1   ! vector opt. (forced unrolling)
402#  else
403!CDIR NOVERRCHK
404      DO jj = 2, jpjm1
405!CDIR NOVERRCHK
406         DO ji = 2, jpim1
407#  endif
408             ikbu = mbku(ji,jj)       ! deepest ocean level at u- and v-points
409             ikbv = mbkv(ji,jj)
410             zfru = 0.5 * fse3u(ji,jj,ikbu) / rdt
411             zfrv = 0.5 * fse3v(ji,jj,ikbv) / rdt
412             IF( ABS( bfrcoef2d(ji,jj) ) > zfru ) THEN
413                IF( ln_ctl ) THEN
414                   WRITE(numout,*) 'BFR ', narea, nimpp+ji, njmpp+jj, ikbu
415                   WRITE(numout,*) 'BFR ', ABS( bfrcoef2d(ji,jj) ), zfru
416                ENDIF
417                ictu = ictu + 1
418             ENDIF
419             IF( ABS( bfrcoef2d(ji,jj) ) > zfrv ) THEN
420                 IF( ln_ctl ) THEN
421                     WRITE(numout,*) 'BFR ', narea, nimpp+ji, njmpp+jj, ikbv
422                     WRITE(numout,*) 'BFR ', bfrcoef2d(ji,jj), zfrv
423                 ENDIF
424                 ictv = ictv + 1
425             ENDIF
426             zminbfr = MIN(  zminbfr, MIN( zfru, ABS( bfrcoef2d(ji,jj) ) )  )
427             zmaxbfr = MAX(  zmaxbfr, MIN( zfrv, ABS( bfrcoef2d(ji,jj) ) )  )
428! (ISF)
429             ikbu = miku(ji,jj)       ! deepest ocean level at u- and v-points
430             ikbv = mikv(ji,jj)
431             zfru = 0.5 * fse3u(ji,jj,ikbu) / rdt
432             zfrv = 0.5 * fse3v(ji,jj,ikbv) / rdt
433             IF( ABS( tfrcoef2d(ji,jj) ) > zfru ) THEN
434                IF( ln_ctl ) THEN
435                   WRITE(numout,*) 'BFR ', narea, nimpp+ji, njmpp+jj, ikbu
436                   WRITE(numout,*) 'BFR ', ABS( tfrcoef2d(ji,jj) ), zfru
437                ENDIF
438                ictu = ictu + 1
439             ENDIF
440             IF( ABS( tfrcoef2d(ji,jj) ) > zfrv ) THEN
441                 IF( ln_ctl ) THEN
442                     WRITE(numout,*) 'BFR ', narea, nimpp+ji, njmpp+jj, ikbv
443                     WRITE(numout,*) 'BFR ', tfrcoef2d(ji,jj), zfrv
444                 ENDIF
445                 ictv = ictv + 1
446             ENDIF
447             zmintfr = MIN(  zmintfr, MIN( zfru, ABS( tfrcoef2d(ji,jj) ) )  )
448             zmaxtfr = MAX(  zmaxtfr, MIN( zfrv, ABS( tfrcoef2d(ji,jj) ) )  )
449
450         END DO
451      END DO
452      IF( lk_mpp ) THEN
453         CALL mpp_sum( ictu )
454         CALL mpp_sum( ictv )
455         CALL mpp_min( zminbfr )
456         CALL mpp_max( zmaxbfr )
457         CALL mpp_min( zmintfr )
458         CALL mpp_max( zmaxtfr )
459      ENDIF
460      IF( .NOT.ln_bfrimp) THEN
461      IF( lwp .AND. ictu + ictv > 0 ) THEN
462         WRITE(numout,*) ' Bottom friction stability check failed at ', ictu, ' U-points '
463         WRITE(numout,*) ' Bottom friction stability check failed at ', ictv, ' V-points '
464         WRITE(numout,*) ' Bottom friction coefficient now ranges from: ', zminbfr, ' to ', zmaxbfr
465         WRITE(numout,*) ' Bottom friction coefficient now ranges from: ', zmintfr, ' to ', zmaxtfr
466         WRITE(numout,*) ' Bottom friction coefficient will be reduced where necessary'
467      ENDIF
468      ENDIF
469      !
470      IF( nn_timing == 1 )  CALL timing_stop('zdf_bfr_init')
471      !
472   END SUBROUTINE zdf_bfr_init
473
474   !!======================================================================
475END MODULE zdfbfr
Note: See TracBrowser for help on using the repository browser.