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/trunk/src/OCE/ISF – NEMO

source: NEMO/trunk/src/OCE/ISF/isfcpl.F90 @ 14143

Last change on this file since 14143 was 14143, checked in by techene, 3 years ago

#2385 add key_linssh equivalent to ln_linssh using domzr_substitute

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