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.
limhdf.F90 in branches/2017/dev_r7881_HPC09_ZDF/NEMOGCM/NEMO/LIM_SRC_3 – NEMO

source: branches/2017/dev_r7881_HPC09_ZDF/NEMOGCM/NEMO/LIM_SRC_3/limhdf.F90 @ 7953

Last change on this file since 7953 was 7753, checked in by mocavero, 7 years ago

Reverting trunk to remove OpenMP

  • Property svn:keywords set to Id
File size: 14.4 KB
Line 
1MODULE limhdf
2   !!======================================================================
3   !!                    ***  MODULE limhdf   ***
4   !! LIM ice model : horizontal diffusion of sea-ice quantities
5   !!======================================================================
6   !! History :  LIM  !  2000-01 (LIM) Original code
7   !!             -   !  2001-05 (G. Madec, R. Hordoir) opa norm
8   !!            1.0  !  2002-08 (C. Ethe)  F90, free form
9   !!            3.6  !  2015-08 (O. Tintó and M. Castrillo)  added lim_hdf (multiple)
10   !!----------------------------------------------------------------------
11#if defined key_lim3
12   !!----------------------------------------------------------------------
13   !!   'key_lim3'                                      LIM3 sea-ice model
14   !!----------------------------------------------------------------------
15   !!   lim_hdf       : diffusion trend on sea-ice variable
16   !!   lim_hdf_init  : initialisation of diffusion trend on sea-ice variable
17   !!----------------------------------------------------------------------
18   USE dom_oce        ! ocean domain
19   USE ice            ! LIM-3: ice variables
20   USE lbclnk         ! lateral boundary condition - MPP exchanges
21   USE lib_mpp        ! MPP library
22   USE wrk_nemo       ! work arrays
23   USE prtctl         ! Print control
24   USE in_out_manager ! I/O manager
25   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
26
27   IMPLICIT NONE
28   PRIVATE
29
30   PUBLIC   lim_hdf         ! called by lim_trp
31   PUBLIC   lim_hdf_init    ! called by sbc_lim_init
32
33   LOGICAL  ::   linit = .TRUE.                             ! initialization flag (set to flase after the 1st call)
34   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   efact   ! metric coefficient
35
36   !! * Substitution
37#  include "vectopt_loop_substitute.h90"
38   !!----------------------------------------------------------------------
39   !! NEMO/LIM3 4.0 , UCL - NEMO Consortium (2010)
40   !! $Id$
41   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
42   !!----------------------------------------------------------------------
43CONTAINS
44
45   SUBROUTINE lim_hdf( ptab, ihdf_vars )
46      !!-------------------------------------------------------------------
47      !!                  ***  ROUTINE lim_hdf  ***
48      !!
49      !! ** purpose :   Compute and add the diffusive trend on sea-ice variables
50      !!
51      !! ** method  :   Second order diffusive operator evaluated using a
52      !!              Cranck-Nicholson time Scheme.
53      !!
54      !! ** Action  :    update ptab with the diffusive contribution
55      !!-------------------------------------------------------------------
56      INTEGER,                    INTENT( in )            ::  ihdf_vars ! number of fields to diffuse
57      REAL(wp), DIMENSION(:,:,:), INTENT( inout ), TARGET ::  ptab      ! Field on which the diffusion is applied
58      !
59      INTEGER                             ::  ji, jj, jk, jl, jm        ! dummy loop indices
60      INTEGER                             ::  iter, ierr, isize         ! local integers
61      REAL(wp)                            ::  zrlxint
62      CHARACTER(lc)                       ::  charout                   ! local character
63      REAL(wp), PARAMETER                 ::  zrelax = 0.5_wp           ! relaxation constant for iterative procedure
64      REAL(wp), PARAMETER                 ::  zalfa  = 0.5_wp           ! =1.0/0.5/0.0 = implicit/Cranck-Nicholson/explicit
65      INTEGER , PARAMETER                 ::  num_iter_max = 100        ! Maximum number of iteration
66      INTEGER , PARAMETER                 ::  num_convfrq  = 5          ! convergence check frequency of the Crant-Nicholson scheme (perf. optimization)
67      REAL(wp), POINTER, DIMENSION(:)     ::  zconv
68      REAL(wp), POINTER, DIMENSION(:,:,:) ::  zrlx, zdiv0, ztab0
69      REAL(wp), POINTER, DIMENSION(:,:)   ::  zflu, zflv, zdiv
70      !!-------------------------------------------------------------------
71      TYPE(arrayptr)   , ALLOCATABLE, DIMENSION(:) ::   pt2d_array, zrlx_array
72      CHARACTER(len=1) , ALLOCATABLE, DIMENSION(:) ::   type_array      ! define the nature of ptab array grid-points
73      !                                                                 ! = T , U , V , F , W and I points
74      REAL(wp)         , ALLOCATABLE, DIMENSION(:) ::   psgn_array      ! =-1 the sign change across the north fold boundary
75      !!-------------------------------------------------------------------
76     
77      !                       !==  Initialisation  ==!
78      ! +1 open water diffusion
79      isize = jpl * ( ihdf_vars + nlay_i ) + 1
80      ALLOCATE( zconv (isize) )
81      ALLOCATE( pt2d_array(isize) , zrlx_array(isize) )
82      ALLOCATE( type_array(isize) )
83      ALLOCATE( psgn_array(isize) )
84
85      CALL wrk_alloc( jpi,jpj,       zflu, zflv, zdiv )
86      CALL wrk_alloc( jpi,jpj,isize, zrlx, zdiv0, ztab0 )
87     
88      DO jk= 1, isize
89         pt2d_array(jk)%pt2d => ptab(:,:,jk)
90         zrlx_array(jk)%pt2d => zrlx(:,:,jk)
91         type_array(jk) = 'T'
92         psgn_array(jk) = 1.
93      END DO
94
95      !
96      IF( linit ) THEN              ! Metric coefficient (compute at the first call and saved in efact)
97         ALLOCATE( efact(jpi,jpj) , STAT=ierr )
98         IF( lk_mpp    )   CALL mpp_sum( ierr )
99         IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'lim_hdf : unable to allocate arrays' )
100         DO jj = 2, jpjm1 
101            DO ji = fs_2 , fs_jpim1   ! vector opt.
102               efact(ji,jj) = ( e2u(ji,jj) + e2u(ji-1,jj) + e1v(ji,jj) + e1v(ji,jj-1) ) * r1_e1e2t(ji,jj)
103            END DO
104         END DO
105         linit = .FALSE.
106      ENDIF
107      !
108      ! Arrays initialization
109      zflu(jpi,:) = 0._wp   
110      zflv(jpi,:) = 0._wp
111      DO jk=1 , isize
112         ztab0(:, : , jk ) = ptab(:,:,jk)
113         zdiv0(:, 1 , jk ) = 0._wp
114         zdiv0(:,jpj, jk ) = 0._wp
115         zdiv0(1,  :, jk ) = 0._wp
116         zdiv0(jpi,:, jk ) = 0._wp
117      END DO
118
119      !                !==  horizontal diffusion using a Crant-Nicholson scheme  ==!
120      zconv(:) = 1._wp
121      iter     = 0
122      !
123      DO WHILE( MAXVAL( zconv(:) ) > ( 2._wp * 1.e-04 ) .AND. iter <= num_iter_max )   ! Sub-time step loop
124         !
125         iter = iter + 1                                 ! incrementation of the sub-time step number
126         !
127         DO jk = 1 , isize
128            jl = ( jk - 1 ) / ( ihdf_vars + nlay_i ) + 1
129            IF ( zconv(jk) > ( 2._wp * 1.e-04 ) ) THEN
130               DO jj = 1, jpjm1                                ! diffusive fluxes in U- and V- direction
131                  DO ji = 1 , fs_jpim1   ! vector opt.
132                     zflu(ji,jj) = pahu3D(ji,jj,jl) * e2u(ji,jj) * r1_e1u(ji,jj) * ( ptab(ji+1,jj,jk) - ptab(ji,jj,jk) )
133                     zflv(ji,jj) = pahv3D(ji,jj,jl) * e1v(ji,jj) * r1_e2v(ji,jj) * ( ptab(ji,jj+1,jk) - ptab(ji,jj,jk) )
134                  END DO
135               END DO
136               !
137               DO jj= 2, jpjm1                                 ! diffusive trend : divergence of the fluxes
138                  DO ji = fs_2 , fs_jpim1   ! vector opt.
139                     zdiv(ji,jj) = ( zflu(ji,jj) - zflu(ji-1,jj) + zflv(ji,jj) - zflv(ji,jj-1) ) * r1_e1e2t(ji,jj)
140                  END DO
141               END DO
142               !
143               IF( iter == 1 )   zdiv0(:,:,jk) = zdiv(:,:)        ! save the 1st evaluation of the diffusive trend in zdiv0
144               !
145               DO jj = 2, jpjm1                                ! iterative evaluation
146                  DO ji = fs_2 , fs_jpim1   ! vector opt.
147                     zrlxint = (   ztab0(ji,jj,jk)    &
148                        &       +  rdt_ice * (           zalfa   * ( zdiv(ji,jj) + efact(ji,jj) * ptab(ji,jj,jk) )   &
149                        &                      + ( 1.0 - zalfa ) *   zdiv0(ji,jj,jk) )                               &
150                        &      ) / ( 1.0 + zalfa * rdt_ice * efact(ji,jj) )
151                     zrlx(ji,jj,jk) = ptab(ji,jj,jk) + zrelax * ( zrlxint - ptab(ji,jj,jk) )
152                  END DO
153               END DO
154            END IF
155
156         END DO
157
158         CALL lbc_lnk_multi( zrlx_array, type_array , psgn_array , isize ) ! Multiple interchange of all the variables
159         !
160
161         IF ( MOD( iter-1 , num_convfrq ) == 0 )  THEN   ! Convergence test every num_convfrq iterations (perf. optimization )
162            DO jk = 1, isize
163               zconv(jk) = 0._wp                                   ! convergence test
164               DO jj = 2, jpjm1
165                  DO ji = fs_2, fs_jpim1
166                     zconv(jk) = MAX( zconv(jk), ABS( zrlx(ji,jj,jk) - ptab(ji,jj,jk) )  )
167                  END DO
168               END DO
169            END DO
170            IF( lk_mpp ) CALL mpp_max_multiple( zconv , isize )            ! max over the global domain for all the variables
171         ENDIF
172         !
173         DO jk=1,isize
174            ptab(:,:,jk) = zrlx(:,:,jk)
175         END DO
176         !
177      END DO  ! end of sub-time step loop
178
179     ! --- final step --- !
180      DO jk = 1, isize
181         jl = ( jk - 1 ) / ( ihdf_vars + nlay_i ) + 1
182         DO jj = 1, jpjm1                                ! diffusive fluxes in U- and V- direction
183            DO ji = 1 , fs_jpim1   ! vector opt.
184               zflu(ji,jj) = pahu3D(ji,jj,jl) * e2u(ji,jj) * r1_e1u(ji,jj) * ( ptab(ji+1,jj,jk) - ptab(ji,jj,jk) )
185               zflv(ji,jj) = pahv3D(ji,jj,jl) * e1v(ji,jj) * r1_e2v(ji,jj) * ( ptab(ji,jj+1,jk) - ptab(ji,jj,jk) )
186            END DO
187         END DO
188         !
189         DO jj= 2, jpjm1                                 ! diffusive trend : divergence of the fluxes
190            DO ji = fs_2 , fs_jpim1   ! vector opt.
191               zdiv(ji,jj) = ( zflu(ji,jj) - zflu(ji-1,jj) + zflv(ji,jj) - zflv(ji,jj-1) ) * r1_e1e2t(ji,jj)
192               ptab(ji,jj,jk) = ztab0(ji,jj,jk) + 0.5 * ( zdiv(ji,jj) + zdiv0(ji,jj,jk) )
193            END DO
194         END DO
195      END DO
196
197      CALL lbc_lnk_multi( pt2d_array, type_array , psgn_array , isize ) ! Multiple interchange of all the variables
198
199      !
200      IF(ln_ctl)   THEN
201         DO jk = 1 , isize
202            zrlx(:,:,jk) = ptab(:,:,jk) - ztab0(:,:,jk)
203            WRITE(charout,FMT="(' lim_hdf  : zconv =',D23.16, ' iter =',I4,2X)") zconv, iter
204            CALL prt_ctl( tab2d_1=zrlx(:,:,jk), clinfo1=charout )
205         END DO
206      ENDIF
207      !
208      CALL wrk_dealloc( jpi,jpj,       zflu, zflv, zdiv )
209      CALL wrk_dealloc( jpi,jpj,isize, zrlx, zdiv0, ztab0 )
210      !
211      DEALLOCATE( zconv )
212      DEALLOCATE( pt2d_array , zrlx_array )
213      DEALLOCATE( type_array )
214      DEALLOCATE( psgn_array )
215      !
216   END SUBROUTINE lim_hdf
217
218   
219   SUBROUTINE lim_hdf_init
220      !!-------------------------------------------------------------------
221      !!                  ***  ROUTINE lim_hdf_init  ***
222      !!
223      !! ** Purpose : Initialisation of horizontal diffusion of sea-ice
224      !!
225      !! ** Method  : Read the namicehdf namelist
226      !!
227      !! ** input   : Namelist namicehdf
228      !!-------------------------------------------------------------------
229      INTEGER  ::   ios                 ! Local integer output status for namelist read
230      NAMELIST/namicehdf/ nn_ahi0, rn_ahi0_ref
231      INTEGER  ::   ji, jj
232      REAL(wp) ::   za00, zd_max
233      !!-------------------------------------------------------------------
234      !
235      REWIND( numnam_ice_ref )              ! Namelist namicehdf in reference namelist : Ice horizontal diffusion
236      READ  ( numnam_ice_ref, namicehdf, IOSTAT = ios, ERR = 901)
237901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicehdf in reference namelist', lwp )
238
239      REWIND( numnam_ice_cfg )              ! Namelist namicehdf in configuration namelist : Ice horizontal diffusion
240      READ  ( numnam_ice_cfg, namicehdf, IOSTAT = ios, ERR = 902 )
241902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicehdf in configuration namelist', lwp )
242      IF(lwm) WRITE ( numoni, namicehdf )
243      !
244      IF(lwp) THEN                          ! control print
245         WRITE(numout,*)
246         WRITE(numout,*) 'lim_hdf_init : Ice horizontal diffusion'
247         WRITE(numout,*) '~~~~~~~~~~~'
248         WRITE(numout,*) '   horizontal diffusivity calculation                          nn_ahi0      = ', nn_ahi0
249         WRITE(numout,*) '   horizontal diffusivity coeff. (orca2 grid)                  rn_ahi0_ref  = ', rn_ahi0_ref
250      ENDIF
251      !
252      !  Diffusion coefficients
253      SELECT CASE( nn_ahi0 )
254
255      CASE( 0 )
256         ahiu(:,:) = rn_ahi0_ref
257         ahiv(:,:) = rn_ahi0_ref
258
259         IF(lwp) WRITE(numout,*) ''
260         IF(lwp) WRITE(numout,*) '   laplacian operator: ahim constant = rn_ahi0_ref'
261
262      CASE( 1 ) 
263
264         zd_max = MAX( MAXVAL( e1t(:,:) ), MAXVAL( e2t(:,:) ) )
265         IF( lk_mpp )   CALL mpp_max( zd_max )          ! max over the global domain
266         
267         ahiu(:,:) = rn_ahi0_ref * zd_max * 1.e-05_wp   ! 1.e05 = 100km = max grid space at 60deg latitude in orca2
268                                                        !                    (60deg = min latitude for ice cover) 
269         ahiv(:,:) = rn_ahi0_ref * zd_max * 1.e-05_wp
270
271         IF(lwp) WRITE(numout,*) ''
272         IF(lwp) WRITE(numout,*) '   laplacian operator: ahim proportional to max of e1 e2 over the domain (', zd_max, ')'
273         IF(lwp) WRITE(numout,*) '   value for ahim = ', rn_ahi0_ref * zd_max * 1.e-05_wp 
274         
275      CASE( 2 ) 
276
277         zd_max = MAX( MAXVAL( e1t(:,:) ), MAXVAL( e2t(:,:) ) )
278         IF( lk_mpp )   CALL mpp_max( zd_max )   ! max over the global domain
279         
280         za00 = rn_ahi0_ref * 1.e-05_wp          ! 1.e05 = 100km = max grid space at 60deg latitude in orca2
281                                                 !                    (60deg = min latitude for ice cover) 
282         DO jj = 1, jpj
283            DO ji = 1, jpi
284               ahiu(ji,jj) = za00 * MAX( e1t(ji,jj), e2t(ji,jj) ) * umask(ji,jj,1)
285               ahiv(ji,jj) = za00 * MAX( e1f(ji,jj), e2f(ji,jj) ) * vmask(ji,jj,1)
286            END DO
287         END DO
288         !
289         IF(lwp) WRITE(numout,*) ''
290         IF(lwp) WRITE(numout,*) '   laplacian operator: ahim proportional to e1'
291         IF(lwp) WRITE(numout,*) '   maximum grid-spacing = ', zd_max, ' maximum value for ahim = ', za00*zd_max
292         
293      END SELECT
294      !
295   END SUBROUTINE lim_hdf_init
296#else
297   !!----------------------------------------------------------------------
298   !!   Default option          Dummy module           NO LIM sea-ice model
299   !!----------------------------------------------------------------------
300#endif
301
302   !!======================================================================
303END MODULE limhdf
Note: See TracBrowser for help on using the repository browser.