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_no_wrk_alloc/NEMOGCM/NEMO/LIM_SRC_3 – NEMO

source: branches/2017/dev_r7881_no_wrk_alloc/NEMOGCM/NEMO/LIM_SRC_3/limhdf.F90 @ 7910

Last change on this file since 7910 was 7910, checked in by timgraham, 7 years ago

All wrk_alloc removed

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