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.
isfcpl.F90 in NEMO/branches/2019/UKMO_MERGE_2019/src/OCE/ISF – NEMO

source: NEMO/branches/2019/UKMO_MERGE_2019/src/OCE/ISF/isfcpl.F90 @ 12077

Last change on this file since 12077 was 12077, checked in by mathiot, 4 years ago

include ENHANCE-02_ISF_nemo in UKMO merge branch

File size: 36.1 KB
Line 
1MODULE isfcpl
2   !!======================================================================
3   !!                       ***  MODULE  isfcpl  ***
4   !!
5   !! iceshelf coupling module : module managing the coupling between NEMO and an ice sheet model
6   !!
7   !!======================================================================
8   !! History :  4.1  !  2019-07  (P. Mathiot) Original code
9   !!----------------------------------------------------------------------
10
11   !!----------------------------------------------------------------------
12   !!   isfrst : read/write iceshelf variables in/from restart
13   !!----------------------------------------------------------------------
14   USE isf_oce                          ! ice shelf variable
15   USE isfutils, ONLY : debug
16   USE lib_mpp , ONLY: mpp_sum, mpp_max ! mpp routine
17   USE domvvl  , ONLY: dom_vvl_zgr      ! vertical scale factor interpolation
18   !
19   USE oce            ! ocean dynamics and tracers
20   USE in_out_manager ! I/O manager
21   USE iom            ! I/O library
22   !
23   IMPLICIT NONE
24
25   PRIVATE
26
27   PUBLIC isfcpl_rst_write, isfcpl_init                    ! iceshelf restart read and write
28   PUBLIC isfcpl_ssh, isfcpl_tra, isfcpl_vol, isfcpl_cons  ! iceshelf correction for ssh, tra, dyn and conservation
29
30   TYPE isfcons
31      INTEGER :: ii     ! i global
32      INTEGER :: jj     ! j global
33      INTEGER :: kk     ! k level
34      REAL(wp):: dvol   ! volume increment
35      REAL(wp):: dsal   ! salt increment
36      REAL(wp):: dtem   ! heat increment
37      REAL(wp):: lon    ! lon
38      REAL(wp):: lat    ! lat
39      INTEGER :: ngb    ! 0/1 (valid location or not (ie on halo or no neigbourg))
40   END TYPE
41   !
42   !!----------------------------------------------------------------------
43   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
44   !! $Id: sbcisf.F90 10536 2019-01-16 19:21:09Z mathiot $
45   !! Software governed by the CeCILL license (see ./LICENSE)
46   !!----------------------------------------------------------------------
47CONTAINS
48   SUBROUTINE isfcpl_init(Kbb, Kmm, Kaa)
49      !!---------------------------------------------------------------------
50      !!                   ***  ROUTINE iscpl_init  ***
51      !!
52      !! ** Purpose : correct ocean state for new wet cell and horizontal divergence
53      !!              correction for the dynamical adjustement
54      !!
55      !! ** Action : - compute ssh on new wet cell
56      !!             - compute T/S on new wet cell
57      !!             - compute horizontal divergence correction as a volume flux
58      !!             - compute the T/S/vol correction increment to keep trend to 0
59      !!
60      !!---------------------------------------------------------------------
61      INTEGER, INTENT(in) :: Kbb, Kmm, Kaa      ! ocean time level indices
62      !!---------------------------------------------------------------------
63      INTEGER :: id
64      !!----------------------------------------------------------------------
65      !
66      ! start on an euler time step
67      neuler = 0
68      !
69      ! allocation and initialisation to 0
70      CALL isf_alloc_cpl()
71      !
72      ! check presence of variable needed for coupling
73      ! iom_varid return 0 if not found
74      id = 1
75      id = id * iom_varid(numror, 'ssmask', ldstop = .false.)
76      id = id * iom_varid(numror, 'tmask' , ldstop = .false.)
77      id = id * iom_varid(numror, 'e3t_n' , ldstop = .false.)
78      id = id * iom_varid(numror, 'e3u_n' , ldstop = .false.)
79      id = id * iom_varid(numror, 'e3v_n' , ldstop = .false.)
80      IF(lwp) WRITE(numout,*) ' isfcpl_init:', id
81      IF (id == 0) THEN
82         IF(lwp) WRITE(numout,*) ' isfcpl_init: restart variables for ice sheet coupling are missing, skip coupling for this leg ' 
83         IF(lwp) WRITE(numout,*) ' ~~~~~~~~~~~'
84         IF(lwp) WRITE(numout,*) ''
85      ELSE
86         ! extrapolation ssh
87         CALL isfcpl_ssh(Kbb, Kmm, Kaa)
88         !
89         ! extrapolation tracer properties
90         CALL isfcpl_tra(Kmm)
91         !
92         ! correction of the horizontal divergence and associated temp. and salt content flux
93         ! Need to : - include in the cpl cons the risfcpl_vol/tsc contribution
94         !           - decide how to manage thickness level change in conservation
95         CALL isfcpl_vol(Kmm)
96         !
97         ! apply the 'conservation' method
98         IF ( ln_isfcpl_cons ) CALL isfcpl_cons(Kmm)
99         !
100      END IF
101      !
102      ! mask velocity properly (mask used in restart not compatible with new mask)
103      uu(:,:,:,Kmm) = uu(:,:,:,Kmm) * umask(:,:,:)
104      vv(:,:,:,Kmm) = vv(:,:,:,Kmm) * vmask(:,:,:)
105      !
106      ! all before fields set to now values
107      ts  (:,:,:,:,Kbb) = ts  (:,:,:,:,Kmm)
108      uu   (:,:,:,Kbb)   = uu   (:,:,:,Kmm)
109      vv   (:,:,:,Kbb)   = vv   (:,:,:,Kmm)
110      ssh (:,:,Kbb)     = ssh (:,:,Kmm)
111      e3t(:,:,:,Kbb)   = e3t(:,:,:,Kmm)
112 
113      ! prepare writing restart
114      IF( lwxios ) THEN
115         CALL iom_set_rstw_var_active('ssmask')
116         CALL iom_set_rstw_var_active('tmask')
117         CALL iom_set_rstw_var_active('e3t_n')
118         CALL iom_set_rstw_var_active('e3u_n')
119         CALL iom_set_rstw_var_active('e3v_n')
120      END IF
121      !
122   END SUBROUTINE isfcpl_init
123   !
124   SUBROUTINE isfcpl_rst_write(kt, Kmm)
125      !!---------------------------------------------------------------------
126      !!                   ***  ROUTINE iscpl_rst_write  ***
127      !!
128      !! ** Purpose : write icesheet coupling variables in restart
129      !!
130      !!-------------------------- IN  --------------------------------------
131      INTEGER, INTENT(in) :: kt
132      INTEGER, INTENT(in) :: Kmm    ! ocean time level index
133      !!----------------------------------------------------------------------
134      !
135      IF( lwxios ) CALL iom_swap( cwxios_context )
136      CALL iom_rstput( kt, nitrst, numrow, 'tmask'  , tmask , ldxios = lwxios )
137      CALL iom_rstput( kt, nitrst, numrow, 'ssmask' , ssmask, ldxios = lwxios )
138      CALL iom_rstput( kt, nitrst, numrow, 'e3t_n'  , e3t(:,:,:,Kmm) , ldxios = lwxios )
139      CALL iom_rstput( kt, nitrst, numrow, 'e3u_n'  , e3u(:,:,:,Kmm) , ldxios = lwxios )
140      CALL iom_rstput( kt, nitrst, numrow, 'e3v_n'  , e3v(:,:,:,Kmm) , ldxios = lwxios )
141      CALL iom_rstput( kt, nitrst, numrow, 'gdepw_n', gdepw(:,:,:,Kmm) , ldxios = lwxios )
142      IF( lwxios ) CALL iom_swap( cxios_context )
143      !
144   END SUBROUTINE isfcpl_rst_write
145
146   SUBROUTINE isfcpl_ssh(Kbb, Kmm, Kaa)
147      !!----------------------------------------------------------------------
148      !!                   ***  ROUTINE iscpl_ssh  ***
149      !!
150      !! ** Purpose :   basic guess of ssh in new wet cell
151      !!
152      !! ** Method  :   basic extrapolation from neigbourg cells
153      !!
154      !!----------------------------------------------------------------------
155      !!
156      INTEGER, INTENT(in) :: Kbb, Kmm, Kaa    ! ocean time level indices
157      !!----------------------------------------------------------------------
158      INTEGER :: ji, jj, jd, jk      !! loop index
159      INTEGER :: jip1, jim1, jjp1, jjm1
160      !!
161      REAL(wp):: zsummsk
162      REAL(wp), DIMENSION(jpi,jpj) :: zdssmask, zssmask0, zssmask_b, zssh
163      !!----------------------------------------------------------------------
164      !
165      CALL iom_get( numror, jpdom_autoglo, 'ssmask'  , zssmask_b, ldxios = lrxios   ) ! need to extrapolate T/S
166
167      ! compute new ssh if we open a full water column
168      ! rude average of the closest neigbourgs (e1e2t not taking into account)
169      !
170      zssh(:,:)     = ssh(:,:,Kmm)
171      zssmask0(:,:) = zssmask_b(:,:)
172      !
173      DO jd = 1, nn_drown
174         !
175         zdssmask(:,:) = ssmask(:,:) - zssmask0(:,:)
176         DO jj = 2,jpj-1
177            DO ji = 2,jpi-1
178               jip1=ji+1; jim1=ji-1;
179               jjp1=jj+1; jjm1=jj-1;
180               !
181               zsummsk = zssmask0(jip1,jj) + zssmask0(jim1,jj) + zssmask0(ji,jjp1) + zssmask0(ji,jjm1)
182               !
183               IF (zdssmask(ji,jj) == 1._wp .AND. zsummsk /= 0._wp) THEN
184                  ssh(ji,jj,Kmm)=( zssh(jip1,jj)*zssmask0(jip1,jj)     &
185                  &           + zssh(jim1,jj)*zssmask0(jim1,jj)     &
186                  &           + zssh(ji,jjp1)*zssmask0(ji,jjp1)     &
187                  &           + zssh(ji,jjm1)*zssmask0(ji,jjm1)) / zsummsk
188                  zssmask_b(ji,jj) = 1._wp
189               ENDIF
190            END DO
191         END DO
192         !
193         zssh(:,:) = ssh(:,:,Kmm)
194         zssmask0(:,:) = zssmask_b(:,:)
195         !
196         CALL lbc_lnk_multi( 'iscplrst', zssh, 'T', 1., zssmask0, 'T', 1. )
197         !
198      END DO
199      !
200      ! update ssh(:,:,Kmm)
201      ssh(:,:,Kmm) = zssh(:,:) * ssmask(:,:)
202      !
203      ssh(:,:,Kbb) = ssh(:,:,Kmm)
204      !
205      IF ( ln_isfdebug ) CALL debug('isfcpl_ssh: sshn',ssh(:,:,Kmm))
206      !
207      ! recompute the vertical scale factor, depth and water thickness
208      IF(lwp) write(numout,*) 'isfcpl_ssh : recompute scale factor from ssh (new wet cell,Kmm)'
209      IF(lwp) write(numout,*) '~~~~~~~~~~~'
210      DO jk = 1, jpk
211         e3t(:,:,jk,Kmm) =  e3t_0(:,:,jk) * ( ht_0(:,:) + ssh(:,:,Kmm) ) &
212             &                          / ( ht_0(:,:) + 1._wp - ssmask(:,:) ) * tmask(:,:,jk)   &
213             &          + e3t_0(:,:,jk)                               * (1._wp -tmask(:,:,jk))
214      END DO
215      e3t(:,:,:,Kbb) = e3t(:,:,:,Kmm)
216      CALL dom_vvl_zgr(Kbb, Kmm, Kaa)
217      !
218   END SUBROUTINE isfcpl_ssh
219
220   SUBROUTINE isfcpl_tra(Kmm)
221      !!----------------------------------------------------------------------
222      !!                   ***  ROUTINE iscpl_tra  ***
223      !!
224      !! ** Purpose :   compute new tn, sn in case of evolving geometry of ice shelves
225      !!
226      !! ** Method  :   tn, sn : basic extrapolation from neigbourg cells
227      !!
228      !!----------------------------------------------------------------------
229      INTEGER, INTENT(in) :: Kmm    ! ocean time level index
230      !!----------------------------------------------------------------------
231      REAL(wp), DIMENSION(jpi,jpj,jpk) :: ztmask_b
232      !REAL(wp), DIMENSION(:,:,:  ), INTENT(in ) :: pdepw_b                         !! depth w before
233      !!
234      INTEGER :: ji, jj, jk, jd          !! loop index
235      INTEGER :: jip1, jim1, jjp1, jjm1, jkp1, jkm1
236      !!
237      REAL(wp):: zsummsk
238      REAL(wp):: zdz, zdzm1, zdzp1
239      !!
240      REAL(wp), DIMENSION(jpi,jpj)          :: zdmask 
241      REAL(wp), DIMENSION(jpi,jpj,jpk)      :: ztmask0, zwmaskn
242      REAL(wp), DIMENSION(jpi,jpj,jpk)      :: ztmask1, zwmaskb, ztmp3d
243      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts) :: zts0
244      !!----------------------------------------------------------------------
245      !
246      CALL iom_get( numror, jpdom_autoglo, 'tmask'  , ztmask_b, ldxios = lrxios   ) ! need to extrapolate T/S
247      !CALL iom_get( numror, jpdom_autoglo, 'wmask'  , zwmask_b, ldxios = lrxios   ) ! need to extrapolate T/S
248      !CALL iom_get( numror, jpdom_autoglo, 'gdepw_n', zdepw_b(:,:,:), ldxios = lrxios ) ! need to interpol vertical profile (vvl)
249      !
250      !
251      ! compute new T/S (interpolation) if vvl only for common wet cell in before and after wmask
252      !PM: Is this IF needed since change to VVL by default
253      !bugged : to be corrected (PM)
254      ! back up original t/s/mask
255      !tsb (:,:,:,:) = ts(:,:,:,:,Kmm)
256      !
257     ! compute new T/S (interpolation) if vvl only for common wet cell in before and after wmask
258
259!      IF (.NOT.ln_linssh) THEN
260!         DO jk = 2,jpk-1
261!            DO jj = 1,jpj
262!               DO ji = 1,jpi
263!                  IF (wmask(ji,jj,jk) * zwmaskb(ji,jj,jk) == 1._wp .AND. (tmask(ji,jj,1)==0._wp .OR. ztmask_b(ji,jj,1)==0._wp) ) THEN
264!
265!                     !compute weight
266!                     zdzp1 = MAX(0._wp,pdepw_b(ji,jj,jk+1) - gdepw(ji,jj,jk+1,Kmm))
267!                     zdzm1 = MAX(0._wp,gdepw(ji,jj,jk  ,Kmm) - pdepw_b(ji,jj,jk  ))
268!                     zdz   = e3t(ji,jj,jk,Kmm) - zdzp1 - zdzm1 ! if isf : e3t = gdepw(ji,jj,jk+1,Kmm)- gdepw(ji,jj,jk,Kmm)
269!
270!                     IF (zdz .LT. 0._wp) THEN
271!                        CALL ctl_stop( 'STOP', 'rst_iscpl : unable to compute the interpolation' )
272!                     END IF
273!
274!                     ts(ji,jj,jk,jp_tem,Kmm) = ( zdzp1*ts(ji,jj,jk+1,jp_tem,Kbb) &
275!                        &                   + zdz  *ts(ji,jj,jk  ,jp_tem,Kbb) &
276!                        &                   + zdzm1*ts(ji,jj,jk-1,jp_tem,Kbb) )/e3t(ji,jj,jk,Kmm)
277!
278!                     ts(ji,jj,jk,jp_sal,Kmm) = ( zdzp1*ts(ji,jj,jk+1,jp_sal,Kbb) &
279!                        &                   + zdz  *ts(ji,jj,jk  ,jp_sal,Kbb) &
280!                        &                   + zdzm1*ts(ji,jj,jk-1,jp_sal,Kbb) )/e3t(ji,jj,jk,Kmm)
281!
282!                  END IF
283!               END DO
284!            END DO
285!         END DO
286!      END IF
287
288      zts0(:,:,:,:)  = ts(:,:,:,:,Kmm)
289      ztmask0(:,:,:) = ztmask_b(:,:,:)
290      ztmask1(:,:,:) = ztmask_b(:,:,:)
291      !
292      ! iterate the extrapolation processes nn_drown times
293      DO jd = 1,nn_drown ! resolution dependent (OK for ISOMIP+ case)
294         DO jk = 1,jpk-1
295            !
296            ! define new wet cell
297            zdmask(:,:) = tmask(:,:,jk) - ztmask0(:,:,jk);
298            !
299            DO jj = 2,jpj-1
300               DO ji = 2,jpi-1
301                  jip1=ji+1; jim1=ji-1;
302                  jjp1=jj+1; jjm1=jj-1;
303                  !
304                  ! check if a wet neigbourg cell is present
305                  zsummsk = ztmask0(jip1,jj  ,jk) + ztmask0(jim1,jj  ,jk) &
306                          + ztmask0(ji  ,jjp1,jk) + ztmask0(ji  ,jjm1,jk)
307                  !
308                  ! if neigbourg wet cell available at the same level
309                  IF ( zdmask(ji,jj) == 1._wp  .AND. zsummsk /= 0._wp ) THEN
310                     !
311                     ! horizontal basic extrapolation
312                     ts(ji,jj,jk,1,Kmm)=( zts0(jip1,jj  ,jk,1) * ztmask0(jip1,jj  ,jk) &
313                     &               + zts0(jim1,jj  ,jk,1) * ztmask0(jim1,jj  ,jk) &
314                     &               + zts0(ji  ,jjp1,jk,1) * ztmask0(ji  ,jjp1,jk) &
315                     &               + zts0(ji  ,jjm1,jk,1) * ztmask0(ji  ,jjm1,jk) ) / zsummsk
316                     ts(ji,jj,jk,2,Kmm)=( zts0(jip1,jj  ,jk,2) * ztmask0(jip1,jj  ,jk) &
317                     &               + zts0(jim1,jj  ,jk,2) * ztmask0(jim1,jj  ,jk) &
318                     &               + zts0(ji  ,jjp1,jk,2) * ztmask0(ji  ,jjp1,jk) &
319                     &               + zts0(ji  ,jjm1,jk,2) * ztmask0(ji  ,jjm1,jk) ) / zsummsk
320                     !
321                     ! update mask for next pass
322                     ztmask1(ji,jj,jk)=1
323                     !
324                  ! in case no neigbourg wet cell available at the same level
325                  ! check if a wet cell is available below
326                  ELSEIF (zdmask(ji,jj) == 1._wp .AND. zsummsk == 0._wp) THEN
327                     !
328                     ! vertical extrapolation if horizontal extrapolation failed
329                     jkm1=max(1,jk-1) ; jkp1=min(jpk,jk+1)
330                     !
331                     ! check if a wet neigbourg cell is present
332                     zsummsk = ztmask0(ji,jj,jkm1) + ztmask0(ji,jj,jkp1)
333                     IF (zdmask(ji,jj) == 1._wp .AND. zsummsk /= 0._wp ) THEN
334                        ts(ji,jj,jk,1,Kmm)=( zts0(ji,jj,jkp1,1)*ztmask0(ji,jj,jkp1)     &
335                        &               + zts0(ji,jj,jkm1,1)*ztmask0(ji,jj,jkm1)) / zsummsk
336                        ts(ji,jj,jk,2,Kmm)=( zts0(ji,jj,jkp1,2)*ztmask0(ji,jj,jkp1)     &
337                        &               + zts0(ji,jj,jkm1,2)*ztmask0(ji,jj,jkm1)) / zsummsk
338                        !
339                        ! update mask for next pass
340                        ztmask1(ji,jj,jk)=1._wp
341                     END IF
342                  END IF
343               END DO
344            END DO
345         END DO
346         !
347         ! update temperature and salinity and mask
348         zts0(:,:,:,:)  = ts(:,:,:,:,Kmm)
349         ztmask0(:,:,:) = ztmask1(:,:,:)
350         !
351         CALL lbc_lnk_multi( 'iscplrst', zts0(:,:,:,jp_tem), 'T', 1., zts0(:,:,:,jp_sal), 'T', 1., ztmask0, 'T', 1.)
352         !
353      END DO  ! nn_drown
354      !
355      ! mask new ts(:,:,:,:,Kmm) field
356      ts(:,:,:,jp_tem,Kmm) = zts0(:,:,:,jp_tem) * tmask(:,:,:)
357      ts(:,:,:,jp_sal,Kmm) = zts0(:,:,:,jp_sal) * tmask(:,:,:)
358      !
359      ! sanity check
360      ! -----------------------------------------------------------------------------------------
361      ! case we open a cell but no neigbour cells available to get an estimate of T and S
362      DO jk = 1,jpk-1
363         DO jj = 1,jpj
364            DO ji = 1,jpi
365               IF (tmask(ji,jj,jk) == 1._wp .AND. ts(ji,jj,jk,2,Kmm) == 0._wp)              &
366                  &   CALL ctl_stop('STOP', 'failing to fill all new weet cell,     &
367                  &                          try increase nn_drown or activate XXXX &
368                  &                         in your domain cfg computation'         )
369            END DO
370         END DO
371      END DO
372      !
373   END SUBROUTINE isfcpl_tra
374
375   SUBROUTINE isfcpl_vol(Kmm)
376      !!----------------------------------------------------------------------
377      !!                   ***  ROUTINE iscpl_vol  ***
378      !!
379      !! ** Purpose : compute the correction of the local divergence to apply 
380      !!              during the first time step after the coupling.
381      !!
382      !! ** Method  : - compute horizontal vol div. before/after coupling
383      !!              - compute vertical input
384      !!              - compute correction
385      !!               
386      !!----------------------------------------------------------------------
387      !!
388      INTEGER, INTENT(in) :: Kmm    ! ocean time level index
389      !!----------------------------------------------------------------------
390      INTEGER :: ji, jj, jk 
391      INTEGER :: ikb, ikt
392      !!
393      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zqvolb, zqvoln  ! vol flux div.         before/after coupling
394      REAL(wp), DIMENSION(jpi,jpj,jpk) :: ze3u_b, ze3v_b  ! vertical scale factor before/after coupling
395      REAL(wp), DIMENSION(jpi,jpj,jpk) :: ztmask_b        ! mask                  before       coupling
396      !!----------------------------------------------------------------------
397      !
398      CALL iom_get( numror, jpdom_autoglo, 'tmask'  , ztmask_b, ldxios = lrxios )
399      CALL iom_get( numror, jpdom_autoglo, 'e3u_n'  , ze3u_b  , ldxios = lrxios )
400      CALL iom_get( numror, jpdom_autoglo, 'e3v_n'  , ze3v_b  , ldxios = lrxios )
401      !
402      ! 1.0: compute horizontal volume flux divergence difference before-after coupling
403      !
404      DO jk = 1, jpk                                 ! Horizontal slab
405         ! 1.1: get volume flux before coupling (>0 out)
406         DO jj = 2, jpjm1
407            DO ji = 2, jpim1
408               zqvolb(ji,jj,jk) =  (   e2u(ji,jj) * ze3u_b(ji,jj,jk) * uu(ji,jj,jk,Kmm) - e2u(ji-1,jj  ) * ze3u_b(ji-1,jj  ,jk) * uu(ji-1,jj  ,jk,Kmm)    &
409                  &                  + e1v(ji,jj) * ze3v_b(ji,jj,jk) * vv(ji,jj,jk,Kmm) - e1v(ji  ,jj-1) * ze3v_b(ji  ,jj-1,jk) * vv(ji  ,jj-1,jk,Kmm)  ) &
410                  &                * ztmask_b(ji,jj,jk)
411            END DO
412         ENDDO
413         !
414         ! 1.2: get volume flux after coupling (>0 out)
415         ! properly mask velocity
416         ! (velocity are still mask with old mask at this stage)
417         uu(:,:,jk,Kmm) = uu(:,:,jk,Kmm) * umask(:,:,jk)
418         vv(:,:,jk,Kmm) = vv(:,:,jk,Kmm) * vmask(:,:,jk)
419         ! compute volume flux divergence after coupling
420         DO jj = 2, jpjm1
421            DO ji = 2, jpim1
422               zqvoln(ji,jj,jk) = (   e2u(ji,jj) * e3u(ji,jj,jk,Kmm) * uu(ji,jj,jk,Kmm) - e2u(ji-1,jj  ) * e3u(ji-1,jj  ,jk,Kmm) * uu(ji-1,jj  ,jk,Kmm)    &
423                  &                 + e1v(ji,jj) * e3v(ji,jj,jk,Kmm) * vv(ji,jj,jk,Kmm) - e1v(ji  ,jj-1) * e3v(ji  ,jj-1,jk,Kmm) * vv(ji  ,jj-1,jk,Kmm)  ) &
424                  &               * tmask(ji,jj,jk)
425            END DO
426         ENDDO
427         !
428         ! 1.3: get 3d volume flux difference (before - after cpl) (>0 out)
429         !      correction to add is _b - _n
430         risfcpl_vol(:,:,jk) = zqvolb(:,:,jk) - zqvoln(:,:,jk)
431      END DO
432      !
433      ! 2.0: include the contribution of the vertical velocity in the volume flux correction
434      !
435      DO jj = 2, jpjm1
436         DO ji = 2, jpim1
437            !
438            ikt = mikt(ji,jj)
439            IF ( ikt > 1 .AND. ssmask(ji,jj) == 1 ) THEN
440               risfcpl_vol(ji,jj,ikt) = risfcpl_vol(ji,jj,ikt) + SUM(zqvolb(ji,jj,1:ikt-1))  ! test sign
441            ENDIF
442            !
443         END DO
444      ENDDO
445      !
446      CALL lbc_lnk( 'iscpl', risfcpl_vol, 'T', 1. )
447      !
448      ! 3.0: set total correction (div, tr(:,:,:,:,Krhs), ssh)
449      !
450      ! 3.1: mask volume flux divergence correction
451      risfcpl_vol(:,:,:) = risfcpl_vol(:,:,:) * tmask(:,:,:)
452      !
453      ! 3.2: get 3d tr(:,:,:,:,Krhs) increment to apply at the first time step
454      ! temperature and salt content flux computed using local ts(:,:,:,:,Kmm)
455      ! (very simple advection scheme)
456      ! (>0 out)
457      risfcpl_tsc(:,:,:,jp_tem) = -risfcpl_vol(:,:,:) * ts(:,:,:,jp_tem,Kmm)
458      risfcpl_tsc(:,:,:,jp_sal) = -risfcpl_vol(:,:,:) * ts(:,:,:,jp_sal,Kmm)
459      !
460      ! 3.3: ssh correction (for dynspg_ts)
461      risfcpl_ssh(:,:) = 0.0
462      DO jk = 1,jpk
463         risfcpl_ssh(:,:) = risfcpl_ssh(:,:) + risfcpl_vol(:,:,jk) * r1_e1e2t(:,:)
464      END DO
465
466   END SUBROUTINE isfcpl_vol
467
468   SUBROUTINE isfcpl_cons(Kmm)
469      !!----------------------------------------------------------------------
470      !!                   ***  ROUTINE iscpl_cons  ***
471      !!
472      !! ** Purpose :   compute the corrective increment in volume/salt/heat to put back the vol/heat/salt
473      !!                removed or added during the coupling processes (wet or dry new cell)
474      !!
475      !! ** Method  :   - compare volume/heat/salt before and after
476      !!                - look for the closest wet cells (share amoung neigbourgs if there are)
477      !!                - build the correction increment to applied at each time step
478      !!               
479      !!----------------------------------------------------------------------
480      !
481      TYPE(isfcons), DIMENSION(:),ALLOCATABLE :: zisfpts ! list of point receiving a correction
482      !
483      !!----------------------------------------------------------------------
484      INTEGER, INTENT(in) :: Kmm    ! ocean time level index
485      !!----------------------------------------------------------------------
486      INTEGER  ::   ji   , jj  , jk  , jproc          ! loop index
487      INTEGER  ::   jip1 , jim1, jjp1, jjm1           ! dummy indices
488      INTEGER  ::   iig  , ijg, ik                    ! dummy indices
489      INTEGER  ::   jisf                              ! start, end and current position in the increment array
490      INTEGER  ::   ingb, ifind                       ! 0/1 target found or need to be found
491      INTEGER  ::   nisfl_area                        ! global number of cell concerned by the wet->dry case
492      INTEGER, DIMENSION(jpnij) :: nisfl              ! local  number of cell concerned by the wet->dry case
493      !
494      REAL(wp) ::   z1_sum, z1_rdtiscpl
495      REAL(wp) ::   zdtem, zdsal, zdvol, zratio       ! tem, sal, vol increment
496      REAL(wp) ::   zlon , zlat                       ! target location 
497      REAL(wp), DIMENSION(jpi,jpj,jpk) :: ztmask_b    ! mask before
498      REAL(wp), DIMENSION(jpi,jpj,jpk) :: ze3t_b      ! scale factor before
499      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zt_b      ! scale factor before
500      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zs_b      ! scale factor before
501      !!----------------------------------------------------------------------
502
503      !==============================================================================
504      ! 1.0: initialisation
505      !==============================================================================
506
507      ! get restart variable
508      CALL iom_get( numror, jpdom_autoglo, 'tmask'  , ztmask_b(:,:,:), ldxios = lrxios   ) ! need to extrapolate T/S
509      CALL iom_get( numror, jpdom_autoglo, 'e3t_n'  , ze3t_b(:,:,:)  , ldxios = lrxios )
510      CALL iom_get( numror, jpdom_autoglo, 'tn'     , zt_b(:,:,:)    , ldxios = lrxios )
511      CALL iom_get( numror, jpdom_autoglo, 'sn'     , zs_b(:,:,:)    , ldxios = lrxios )
512
513      ! compute run length
514      nstp_iscpl  = nitend - nit000 + 1
515      rdt_iscpl   = nstp_iscpl * rn_rdt
516      z1_rdtiscpl = 1._wp / rdt_iscpl 
517
518      IF (lwp) WRITE(numout,*) '            nb of stp for cons  = ', nstp_iscpl
519      IF (lwp) WRITE(numout,*) '            coupling time step  = ', rdt_iscpl
520
521      ! initialisation correction
522      risfcpl_cons_vol = 0.0
523      risfcpl_cons_ssh = 0.0
524      risfcpl_cons_tsc = 0.0
525
526      !==============================================================================
527      ! 2.0: diagnose the heat, salt and volume input and compute the correction variable
528      !      for case where we wet a cell or cell still wet (no change in cell status)
529      !==============================================================================
530
531      DO jk = 1,jpk-1
532         DO jj = nldj,nlej
533            DO ji = nldi,nlei
534
535               ! volume diff
536               zdvol = e3t(ji,jj,jk,Kmm) * tmask(ji,jj,jk) - ze3t_b(ji,jj,jk) * ztmask_b(ji,jj,jk)
537
538               ! heat diff
539               zdtem = ts (ji,jj,jk,jp_tem,Kmm) *  e3t(ji,jj,jk,Kmm) *  tmask  (ji,jj,jk)   &
540                     - zt_b(ji,jj,jk)        * ze3t_b(ji,jj,jk) * ztmask_b(ji,jj,jk)
541
542               ! salt diff
543               zdsal = ts(ji,jj,jk,jp_sal,Kmm) *  e3t(ji,jj,jk,Kmm) *  tmask  (ji,jj,jk)   &
544                     - zs_b(ji,jj,jk)       * ze3t_b(ji,jj,jk) * ztmask_b(ji,jj,jk)
545           
546               ! volume, heat and salt differences in each cell (>0 means correction is an outward flux)
547               ! in addition to the geometry change unconservation, need to add the divergence correction as it is flux across the boundary
548               risfcpl_cons_vol(ji,jj,jk)        = (   zdvol * e1e2t(ji,jj) + risfcpl_vol(ji,jj,jk)        ) * z1_rdtiscpl
549               risfcpl_cons_tsc(ji,jj,jk,jp_sal) = ( - zdsal * e1e2t(ji,jj) + risfcpl_tsc(ji,jj,jk,jp_sal) ) * z1_rdtiscpl
550               risfcpl_cons_tsc(ji,jj,jk,jp_tem) = ( - zdtem * e1e2t(ji,jj) + risfcpl_tsc(ji,jj,jk,jp_tem) ) * z1_rdtiscpl
551
552            END DO
553         END DO
554      END DO
555      !
556      !==============================================================================
557      ! 3.0: diagnose the heat, salt and volume input and compute the correction variable
558      !      for case where we close a cell
559      !==============================================================================
560      !
561      ! compute the total number of point receiving a correction increment for each processor
562      ! local
563      nisfl(:)=0
564      DO jk = 1,jpk-1
565         DO jj = nldj,nlej
566            DO ji = nldi,nlei
567               jip1=MIN(ji+1,jpi) ; jim1=MAX(ji-1,1) ; jjp1=MIN(jj+1,jpj) ; jjm1=MAX(jj-1,1) ;
568               IF ( tmask(ji,jj,jk) == 0._wp .AND. ztmask_b(ji,jj,jk) == 1._wp ) nisfl(narea) = nisfl(narea) + MAX(SUM(tmask(jim1:jip1,jjm1:jjp1,jk)),1._wp)
569            ENDDO
570         ENDDO
571      ENDDO
572      !
573      ! global
574      CALL mpp_sum('isfcpl',nisfl  )
575      !
576      ! allocate list of point receiving correction
577      ALLOCATE(zisfpts(nisfl(narea)))
578      !
579      zisfpts(:) = isfcons(0,0,0,-HUGE(1.0), -HUGE(1.0), -HUGE(1.0), -HUGE(1.0), -HUGE(1.0), 0)
580      !
581      ! start computing the correction and fill zisfpts
582      ! local
583      jisf = 0
584      DO jk = 1,jpk-1
585         DO jj = nldj,nlej
586            DO ji = nldi,nlei
587               IF ( tmask(ji,jj,jk) == 0._wp .AND. ztmask_b(ji,jj,jk) == 1._wp ) THEN
588
589                  jip1=MIN(ji+1,jpi) ; jim1=MAX(ji-1,1) ; jjp1=MIN(jj+1,jpj) ; jjm1=MAX(jj-1,1) ;
590
591                  zdvol = risfcpl_cons_vol(ji,jj,jk       )
592                  zdsal = risfcpl_cons_tsc(ji,jj,jk,jp_sal)
593                  zdtem = risfcpl_cons_tsc(ji,jj,jk,jp_tem)
594
595                  IF ( SUM( tmask(jim1:jip1,jjm1:jjp1,jk) ) > 0._wp ) THEN
596                     ! spread correction amoung neigbourg wet cells (horizontal direction first)
597                     ! as it is a rude correction corner and lateral cell have the same weight
598                     !
599                     z1_sum =  1._wp / SUM( tmask(jim1:jip1,jjm1:jjp1,jk) )
600                     !
601                     ! lateral cells
602                     IF (tmask(jip1,jj  ,jk) == 1) CALL update_isfpts(zisfpts, jisf, jip1, jj  , jk, zdvol, zdsal, zdtem, z1_sum)
603                     IF (tmask(jim1,jj  ,jk) == 1) CALL update_isfpts(zisfpts, jisf, jim1, jj  , jk, zdvol, zdsal, zdtem, z1_sum)
604                     IF (tmask(ji  ,jjp1,jk) == 1) CALL update_isfpts(zisfpts, jisf, ji  , jjp1, jk, zdvol, zdsal, zdtem, z1_sum)
605                     IF (tmask(ji  ,jjm1,jk) == 1) CALL update_isfpts(zisfpts, jisf, ji  , jjm1, jk, zdvol, zdsal, zdtem, z1_sum)
606                     !
607                     ! corner  cells
608                     IF (tmask(jip1,jjm1,jk) == 1) CALL update_isfpts(zisfpts, jisf, jip1, jjm1, jk, zdvol, zdsal, zdtem, z1_sum)
609                     IF (tmask(jim1,jjm1,jk) == 1) CALL update_isfpts(zisfpts, jisf, jim1, jjm1, jk, zdvol, zdsal, zdtem, z1_sum)
610                     IF (tmask(jim1,jjp1,jk) == 1) CALL update_isfpts(zisfpts, jisf, jim1, jjp1, jk, zdvol, zdsal, zdtem, z1_sum)
611                     IF (tmask(jip1,jjp1,jk) == 1) CALL update_isfpts(zisfpts, jisf, jip1, jjp1, jk, zdvol, zdsal, zdtem, z1_sum)
612                     !
613                  ELSE IF ( tmask(ji,jj,jk+1) == 1._wp ) THEN
614                     ! spread correction amoung neigbourg wet cells (vertical direction)
615                     CALL update_isfpts(zisfpts, jisf, ji  , jj  , jk+1, zdvol, zdsal, zdtem, 1., 0)
616                  ELSE
617                     ! need to find where to put correction in later on
618                     CALL update_isfpts(zisfpts, jisf, ji  , jj  , jk  , zdvol, zdsal, zdtem, 1., 1)
619                  END IF
620               END IF
621            END DO
622         END DO
623      END DO
624      !
625      ! share data among all processes because for some point we need to find the closest wet point (could be on other process)
626      DO jproc=1,jpnij
627         !
628         ! share total number of isf point treated for proc jproc
629         IF (jproc==narea) THEN
630            nisfl_area=nisfl(jproc)
631         ELSE
632            nisfl_area=0
633         END IF
634         CALL mpp_max('isfcpl',nisfl_area)
635         !
636         DO jisf = 1,nisfl_area
637            !
638            IF (jproc==narea) THEN
639               ! indices (conversion to global indices and sharing)
640               iig = zisfpts(jisf)%ii       ; ijg = zisfpts(jisf)%jj       ; ik = zisfpts(jisf)%kk
641               !
642               ! data
643               zdvol = zisfpts(jisf)%dvol   ; zdsal = zisfpts(jisf)%dsal   ; zdtem = zisfpts(jisf)%dtem
644               !
645               ! location
646               zlat = zisfpts(jisf)%lat     ; zlon = zisfpts(jisf)%lon
647               !
648               ! find flag
649               ingb = zisfpts(jisf)%ngb
650            ELSE
651               iig  =0   ; ijg  =0   ; ik   =0 
652               zdvol=-HUGE(1.0) ; zdsal=-HUGE(1.0) ; zdtem=-HUGE(1.0)
653               zlat =-HUGE(1.0) ; zlon =-HUGE(1.0)   
654               ingb = 0
655            END IF
656            !
657            ! share data (need synchronisation of data as get_correction call a global com)
658            CALL mpp_max('isfcpl',iig)   ; CALL mpp_max('isfcpl',ijg)   ; CALL mpp_max('isfcpl',ik)
659            CALL mpp_max('isfcpl',zdvol) ; CALL mpp_max('isfcpl',zdsal) ; CALL mpp_max('isfcpl',zdtem)
660            CALL mpp_max('isfcpl',zlat)  ; CALL mpp_max('isfcpl',zlon)
661            CALL mpp_max('isfcpl',ingb)
662            !
663            ! fill the 3d correction array
664            CALL get_correction(iig, ijg, ik, zlon, zlat, zdvol, zdsal, zdtem, ingb)
665         END DO
666      END DO
667      !
668      !==============================================================================
669      ! 4.0: finalisation and compute ssh equivalent of the volume correction
670      !==============================================================================
671      !
672      ! mask (>0 out)
673      risfcpl_cons_vol(:,:,:       ) = risfcpl_cons_vol(:,:,:       ) * tmask(:,:,:)
674      risfcpl_cons_tsc(:,:,:,jp_sal) = risfcpl_cons_tsc(:,:,:,jp_sal) * tmask(:,:,:)
675      risfcpl_cons_tsc(:,:,:,jp_tem) = risfcpl_cons_tsc(:,:,:,jp_tem) * tmask(:,:,:)
676      !
677      ! add lbclnk
678      CALL lbc_lnk_multi( 'iscplrst', risfcpl_cons_tsc(:,:,:,jp_tem), 'T', 1., risfcpl_cons_tsc(:,:,:,jp_sal), 'T', 1., &
679         &                            risfcpl_cons_vol(:,:,:)       , 'T', 1.)
680      !
681      ! ssh correction (for dynspg_ts)
682      DO jk = 1,jpk
683         risfcpl_cons_ssh(:,:) = risfcpl_cons_ssh(:,:) + risfcpl_cons_vol(:,:,jk)
684      END DO
685      risfcpl_cons_ssh(:,:) = risfcpl_cons_ssh(:,:) * r1_e1e2t(:,:)
686      !
687   END SUBROUTINE isfcpl_cons
688   !
689   SUBROUTINE update_isfpts(sisfpts, kpts, ki, kj, kk, pdvol, pdsal, pdtem, pratio, kfind)
690      !!---------------------------------------------------------------------
691      !!                  ***  ROUTINE update_isfpts  ***
692      !!
693      !! ** Purpose : if a cell become dry, we need to put the corrective increment elsewhere
694      !!
695      !! ** Action  : update the list of point
696      !!
697      !!----------------------------------------------------------------------
698      !!----------------------------------------------------------------------
699      TYPE(isfcons), DIMENSION(:), INTENT(inout) :: sisfpts
700      INTEGER,                     INTENT(inout) :: kpts
701      !!----------------------------------------------------------------------
702      INTEGER,      INTENT(in   )           :: ki, kj, kk                  !    target location (kfind=0)
703      !                                                                    ! or source location (kfind=1)
704      INTEGER,      INTENT(in   ), OPTIONAL :: kfind                       ! 0  target cell already found
705      !                                                                    ! 1  target to be determined
706      REAL(wp),     INTENT(in   )           :: pdvol, pdsal, pdtem, pratio ! vol/sal/tem increment
707      !                                                                    ! and ratio in case increment span over multiple cells.
708      !!----------------------------------------------------------------------
709      INTEGER :: ifind
710      !!----------------------------------------------------------------------
711      !
712      ! increment position
713      kpts = kpts + 1
714      !
715      ! define if we need to look for closest valid wet cell (no neighbours or neigbourg on halo)
716      IF ( PRESENT(kfind) ) THEN
717         ifind = kfind
718      ELSE
719         ifind = ( 1 - tmask_h(ki,kj) ) * tmask(ki,kj,kk)
720      END IF
721      !
722      ! update isfpts structure
723      sisfpts(kpts) = isfcons(mig(ki), mjg(kj), kk, pratio * pdvol, pratio * pdsal, pratio * pdtem, glamt(ki,kj), gphit(ki,kj), ifind )
724      !
725   END SUBROUTINE update_isfpts
726   !
727   SUBROUTINE get_correction( ki, kj, kk, plon, plat, pvolinc, psalinc, pteminc, kfind)
728      !!---------------------------------------------------------------------
729      !!                  ***  ROUTINE get_correction  ***
730      !!
731      !! ** Action : - Find the closest valid cell if needed (wet and not on the halo)
732      !!             - Scale the correction depending of pratio (case where multiple wet neigbourgs)
733      !!             - Fill the correction array
734      !!
735      !!----------------------------------------------------------------------
736      INTEGER , INTENT(in) :: ki, kj, kk, kfind        ! target point indices
737      REAL(wp), INTENT(in) :: plon, plat               ! target point lon/lat
738      REAL(wp), INTENT(in) :: pvolinc, pteminc,psalinc ! correction increment for vol/temp/salt
739      !!----------------------------------------------------------------------
740      INTEGER :: jj, ji, iig, ijg
741      !!----------------------------------------------------------------------
742      !
743      ! define global indice of correction location
744      iig = ki ; ijg = kj
745      IF ( kfind == 1 ) CALL dom_ngb( plon, plat, iig, ijg,'T', kk)
746      !
747      ! fill the correction array
748      DO jj = mj0(ijg),mj1(ijg)
749         DO ji = mi0(iig),mi1(iig)
750            ! correct the vol_flx and corresponding heat/salt flx in the closest cell
751            risfcpl_cons_vol(ji,jj,kk)        =  risfcpl_cons_vol(ji,jj,kk       ) + pvolinc
752            risfcpl_cons_tsc(ji,jj,kk,jp_sal) =  risfcpl_cons_tsc(ji,jj,kk,jp_sal) + psalinc
753            risfcpl_cons_tsc(ji,jj,kk,jp_tem) =  risfcpl_cons_tsc(ji,jj,kk,jp_tem) + pteminc
754         END DO
755      END DO
756
757   END SUBROUTINE get_correction
758
759END MODULE isfcpl
Note: See TracBrowser for help on using the repository browser.