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.
sbcice_cice.F90 in NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/SBC – NEMO

source: NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/SBC/sbcice_cice.F90 @ 11960

Last change on this file since 11960 was 11960, checked in by acc, 4 years ago

Branch 2019/dev_r11943_MERGE_2019. Merge in changes from 2019/dev_r11613_ENHANCE-04_namelists_as_internalfiles. (svn merge -r 11614:11954). Resolved tree conflicts and one actual conflict. Sette tested(these changes alter the ext/AGRIF reference; remember to update). See ticket #2341

  • Property svn:keywords set to Id
File size: 43.6 KB
Line 
1MODULE sbcice_cice
2   !!======================================================================
3   !!                       ***  MODULE  sbcice_cice  ***
4   !! To couple with sea ice model CICE (LANL)
5   !!=====================================================================
6#if defined key_cice
7   !!----------------------------------------------------------------------
8   !!   'key_cice' :                                     CICE sea-ice model
9   !!----------------------------------------------------------------------
10   !!   sbc_ice_cice  : sea-ice model time-stepping and update ocean sbc over ice-covered area
11   !!----------------------------------------------------------------------
12   USE oce             ! ocean dynamics and tracers
13   USE dom_oce         ! ocean space and time domain
14   USE domvvl
15   USE phycst, only : rcp, rau0, r1_rau0, rhos, rhoi
16   USE in_out_manager  ! I/O manager
17   USE iom, ONLY : iom_put,iom_use              ! I/O manager library !!Joakim edit
18   USE lib_mpp         ! distributed memory computing library
19   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
20   USE daymod          ! calendar
21   USE fldread         ! read input fields
22   USE sbc_oce         ! Surface boundary condition: ocean fields
23   USE sbc_ice         ! Surface boundary condition: ice   fields
24   USE sbcblk          ! Surface boundary condition: bulk
25   USE sbccpl
26
27   USE ice_kinds_mod
28   USE ice_blocks
29   USE ice_domain
30   USE ice_domain_size
31   USE ice_boundary
32   USE ice_constants
33   USE ice_gather_scatter
34   USE ice_calendar, only: dt
35   USE ice_state, only: aice,aicen,uvel,vvel,vsno,vsnon,vice,vicen
36# if defined key_cice4
37   USE ice_flux, only: strax,stray,strocnx,strocny,frain,fsnow,  &
38                strocnxT,strocnyT,                               & 
39                sst,sss,uocn,vocn,ss_tltx,ss_tlty,fsalt_gbm,     &
40                fresh_gbm,fhocn_gbm,fswthru_gbm,frzmlt,          &
41                flatn_f,fsurfn_f,fcondtopn_f,                    &
42                uatm,vatm,wind,fsw,flw,Tair,potT,Qa,rhoa,zlvl,   &
43                swvdr,swvdf,swidr,swidf
44   USE ice_therm_vertical, only: calc_Tsfc
45#else
46   USE ice_flux, only: strax,stray,strocnx,strocny,frain,fsnow,  &
47                strocnxT,strocnyT,                               & 
48                sst,sss,uocn,vocn,ss_tltx,ss_tlty,fsalt_ai,     &
49                fresh_ai,fhocn_ai,fswthru_ai,frzmlt,          &
50                flatn_f,fsurfn_f,fcondtopn_f,                    &
51                uatm,vatm,wind,fsw,flw,Tair,potT,Qa,rhoa,zlvl,   &
52                swvdr,swvdf,swidr,swidf
53   USE ice_therm_shared, only: calc_Tsfc
54#endif
55   USE ice_forcing, only: frcvdr,frcvdf,frcidr,frcidf
56   USE ice_atmo, only: calc_strair
57
58   USE CICE_InitMod
59   USE CICE_RunMod
60   USE CICE_FinalMod
61
62   IMPLICIT NONE
63   PRIVATE
64
65   PUBLIC cice_sbc_init   ! routine called by sbc_init
66   PUBLIC cice_sbc_final  ! routine called by sbc_final
67   PUBLIC sbc_ice_cice    ! routine called by sbc
68
69   INTEGER             ::   ji_off
70   INTEGER             ::   jj_off
71
72   INTEGER , PARAMETER ::   jpfld   = 13   ! maximum number of files to read
73   INTEGER , PARAMETER ::   jp_snow = 1    ! index of snow file
74   INTEGER , PARAMETER ::   jp_rain = 2    ! index of rain file
75   INTEGER , PARAMETER ::   jp_sblm = 3    ! index of sublimation file
76   INTEGER , PARAMETER ::   jp_top1 = 4    ! index of category 1 topmelt file
77   INTEGER , PARAMETER ::   jp_top2 = 5    ! index of category 2 topmelt file
78   INTEGER , PARAMETER ::   jp_top3 = 6    ! index of category 3 topmelt file
79   INTEGER , PARAMETER ::   jp_top4 = 7    ! index of category 4 topmelt file
80   INTEGER , PARAMETER ::   jp_top5 = 8    ! index of category 5 topmelt file
81   INTEGER , PARAMETER ::   jp_bot1 = 9    ! index of category 1 botmelt file
82   INTEGER , PARAMETER ::   jp_bot2 = 10   ! index of category 2 botmelt file
83   INTEGER , PARAMETER ::   jp_bot3 = 11   ! index of category 3 botmelt file
84   INTEGER , PARAMETER ::   jp_bot4 = 12   ! index of category 4 botmelt file
85   INTEGER , PARAMETER ::   jp_bot5 = 13   ! index of category 5 botmelt file
86   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf    ! structure of input fields (file informations, fields read)
87
88   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:), PRIVATE ::   png     ! local array used in sbc_cice_ice
89
90   !!----------------------------------------------------------------------
91   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
92   !! $Id$
93   !! Software governed by the CeCILL license (see ./LICENSE)
94   !!----------------------------------------------------------------------
95CONTAINS
96
97   INTEGER FUNCTION sbc_ice_cice_alloc()
98      !!----------------------------------------------------------------------
99      !!                ***  FUNCTION sbc_ice_cice_alloc  ***
100      !!----------------------------------------------------------------------
101      ALLOCATE( png(jpi,jpj,jpnij), STAT=sbc_ice_cice_alloc )
102      CALL mpp_sum ( 'sbcice_cice', sbc_ice_cice_alloc )
103      IF( sbc_ice_cice_alloc > 0 )   CALL ctl_warn('sbc_ice_cice_alloc: allocation of arrays failed.')
104   END FUNCTION sbc_ice_cice_alloc
105
106   SUBROUTINE sbc_ice_cice( kt, ksbc )
107      !!---------------------------------------------------------------------
108      !!                  ***  ROUTINE sbc_ice_cice  ***
109      !!                   
110      !! ** Purpose :   update the ocean surface boundary condition via the
111      !!                CICE Sea Ice Model time stepping
112      !!
113      !! ** Method  : - Get any extra forcing fields for CICE 
114      !!              - Prepare forcing fields
115      !!              - CICE model time stepping
116      !!              - call the routine that computes mass and
117      !!                heat fluxes at the ice/ocean interface
118      !!
119      !! ** Action  : - time evolution of the CICE sea-ice model
120      !!              - update all sbc variables below sea-ice:
121      !!                utau, vtau, qns , qsr, emp , sfx
122      !!---------------------------------------------------------------------
123      INTEGER, INTENT(in) ::   kt      ! ocean time step
124      INTEGER, INTENT(in) ::   ksbc    ! surface forcing type
125      !!----------------------------------------------------------------------
126      !
127      !                                        !----------------------!
128      IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN     !  Ice time-step only  !
129         !                                     !----------------------!
130
131         ! Make sure any fluxes required for CICE are set
132         IF      ( ksbc == jp_flx ) THEN
133            CALL cice_sbc_force(kt)
134         ELSE IF ( ksbc == jp_purecpl ) THEN
135            CALL sbc_cpl_ice_flx( fr_i )
136         ENDIF
137
138         CALL cice_sbc_in  ( kt, ksbc )
139         CALL CICE_Run
140         CALL cice_sbc_out ( kt, ksbc )
141
142         IF ( ksbc == jp_purecpl )  CALL cice_sbc_hadgam(kt+1)
143
144      ENDIF                                          ! End sea-ice time step only
145      !
146   END SUBROUTINE sbc_ice_cice
147
148
149   SUBROUTINE cice_sbc_init( ksbc, Kbb, Kmm )
150      !!---------------------------------------------------------------------
151      !!                    ***  ROUTINE cice_sbc_init  ***
152      !! ** Purpose: Initialise ice related fields for NEMO and coupling
153      !!
154      !!---------------------------------------------------------------------
155      INTEGER, INTENT( in  ) ::   ksbc                ! surface forcing type
156      INTEGER, INTENT( in  ) ::   Kbb, Kmm            ! time level indices
157      REAL(wp), DIMENSION(jpi,jpj) :: ztmp1, ztmp2
158      REAL(wp) ::   zcoefu, zcoefv, zcoeff            ! local scalar
159      INTEGER  ::   ji, jj, jl, jk                    ! dummy loop indices
160      !!---------------------------------------------------------------------
161      !
162      IF(lwp) WRITE(numout,*)'cice_sbc_init'
163
164      ji_off = INT ( (jpiglo - nx_global) / 2 )
165      jj_off = INT ( (jpjglo - ny_global) / 2 )
166
167#if defined key_nemocice_decomp
168      ! Pass initial SST from NEMO to CICE so ice is initialised correctly if
169      ! there is no restart file.
170      ! Values from a CICE restart file would overwrite this
171      IF ( .NOT. ln_rstart ) THEN   
172         CALL nemo2cice( ts(:,:,1,jp_tem,Kmm) , sst , 'T' , 1.) 
173      ENDIF 
174#endif
175
176! Initialize CICE
177      CALL CICE_Initialize
178
179! Do some CICE consistency checks
180      IF ( (ksbc == jp_flx) .OR. (ksbc == jp_purecpl) ) THEN
181         IF ( calc_strair .OR. calc_Tsfc ) THEN
182            CALL ctl_stop( 'STOP', 'cice_sbc_init : Forcing option requires calc_strair=F and calc_Tsfc=F in ice_in' )
183         ENDIF
184      ELSEIF (ksbc == jp_blk) THEN
185         IF ( .NOT. (calc_strair .AND. calc_Tsfc) ) THEN
186            CALL ctl_stop( 'STOP', 'cice_sbc_init : Forcing option requires calc_strair=T and calc_Tsfc=T in ice_in' )
187         ENDIF
188      ENDIF
189
190
191! allocate sbc_ice and sbc_cice arrays
192      IF( sbc_ice_alloc()      /= 0 )   CALL ctl_stop( 'STOP', 'sbc_ice_cice_alloc : unable to allocate arrays' )
193      IF( sbc_ice_cice_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'sbc_ice_cice_alloc : unable to allocate cice arrays' )
194
195! Ensure ocean temperatures are nowhere below freezing if not a NEMO restart
196      IF( .NOT. ln_rstart ) THEN
197         ts(:,:,:,jp_tem,Kmm) = MAX (ts(:,:,:,jp_tem,Kmm),Tocnfrz)
198         ts(:,:,:,jp_tem,Kbb) = ts(:,:,:,jp_tem,Kmm)
199      ENDIF
200
201      fr_iu(:,:)=0.0
202      fr_iv(:,:)=0.0
203
204      CALL cice2nemo(aice,fr_i, 'T', 1. )
205      IF ( (ksbc == jp_flx) .OR. (ksbc == jp_purecpl) ) THEN
206         DO jl=1,ncat
207            CALL cice2nemo(aicen(:,:,jl,:),a_i(:,:,jl), 'T', 1. )
208         ENDDO
209      ENDIF
210
211! T point to U point
212! T point to V point
213      DO jj=1,jpjm1
214         DO ji=1,jpim1
215            fr_iu(ji,jj)=0.5*(fr_i(ji,jj)+fr_i(ji+1,jj))*umask(ji,jj,1)
216            fr_iv(ji,jj)=0.5*(fr_i(ji,jj)+fr_i(ji,jj+1))*vmask(ji,jj,1)
217         ENDDO
218      ENDDO
219
220      CALL lbc_lnk_multi( 'sbcice_cice', fr_iu , 'U', 1.,  fr_iv , 'V', 1. )
221
222      ! set the snow+ice mass
223      CALL cice2nemo(vsno(:,:,:),ztmp1,'T', 1. )
224      CALL cice2nemo(vice(:,:,:),ztmp2,'T', 1. )
225      snwice_mass  (:,:) = ( rhos * ztmp1(:,:) + rhoi * ztmp2(:,:)  )
226      snwice_mass_b(:,:) = snwice_mass(:,:)
227
228      IF( .NOT.ln_rstart ) THEN
229         IF( ln_ice_embd ) THEN            ! embedded sea-ice: deplete the initial ssh below sea-ice area
230            ssh(:,:,Kmm) = ssh(:,:,Kmm) - snwice_mass(:,:) * r1_rau0
231            ssh(:,:,Kbb) = ssh(:,:,Kbb) - snwice_mass(:,:) * r1_rau0
232
233!!gm This should be put elsewhere....   (same remark for limsbc)
234!!gm especially here it is assumed zstar coordinate, but it can be ztilde....
235            IF( .NOT.ln_linssh ) THEN
236               !
237               DO jk = 1,jpkm1                     ! adjust initial vertical scale factors
238                  e3t(:,:,jk,Kmm) = e3t_0(:,:,jk)*( 1._wp + ssh(:,:,Kmm)*tmask(:,:,1)/(ht_0(:,:) + 1.0 - tmask(:,:,1)) )
239                  e3t(:,:,jk,Kbb) = e3t_0(:,:,jk)*( 1._wp + ssh(:,:,Kbb)*tmask(:,:,1)/(ht_0(:,:) + 1.0 - tmask(:,:,1)) )
240               ENDDO
241               e3t(:,:,:,Krhs) = e3t(:,:,:,Kbb)
242               ! Reconstruction of all vertical scale factors at now and before time-steps
243               ! =============================================================================
244               ! Horizontal scale factor interpolations
245               ! --------------------------------------
246               CALL dom_vvl_interpol( e3t(:,:,:,Kbb), e3u(:,:,:,Kbb), 'U' )
247               CALL dom_vvl_interpol( e3t(:,:,:,Kbb), e3v(:,:,:,Kbb), 'V' )
248               CALL dom_vvl_interpol( e3t(:,:,:,Kmm), e3u(:,:,:,Kmm), 'U' )
249               CALL dom_vvl_interpol( e3t(:,:,:,Kmm), e3v(:,:,:,Kmm), 'V' )
250               CALL dom_vvl_interpol( e3u(:,:,:,Kmm), e3f(:,:,:), 'F' )
251               ! Vertical scale factor interpolations
252               ! ------------------------------------
253               CALL dom_vvl_interpol( e3t(:,:,:,Kmm), e3w (:,:,:,Kmm), 'W'  )
254               CALL dom_vvl_interpol( e3u(:,:,:,Kmm), e3uw(:,:,:,Kmm), 'UW' )
255               CALL dom_vvl_interpol( e3v(:,:,:,Kmm), e3vw(:,:,:,Kmm), 'VW' )
256               CALL dom_vvl_interpol( e3u(:,:,:,Kbb), e3uw(:,:,:,Kbb), 'UW' )
257               CALL dom_vvl_interpol( e3v(:,:,:,Kbb), e3vw(:,:,:,Kbb), 'VW' )
258               ! t- and w- points depth
259               ! ----------------------
260               gdept(:,:,1,Kmm) = 0.5_wp * e3w(:,:,1,Kmm)
261               gdepw(:,:,1,Kmm) = 0.0_wp
262               gde3w(:,:,1)     = gdept(:,:,1,Kmm) - ssh(:,:,Kmm)
263               DO jk = 2, jpk
264                  gdept(:,:,jk,Kmm) = gdept(:,:,jk-1,Kmm) + e3w(:,:,jk,Kmm)
265                  gdepw(:,:,jk,Kmm) = gdepw(:,:,jk-1,Kmm) + e3t(:,:,jk-1,Kmm)
266                  gde3w(:,:,jk)     = gdept(:,:,jk  ,Kmm) - sshn   (:,:)
267               END DO
268            ENDIF
269         ENDIF
270      ENDIF
271      !
272   END SUBROUTINE cice_sbc_init
273
274   
275   SUBROUTINE cice_sbc_in( kt, ksbc )
276      !!---------------------------------------------------------------------
277      !!                    ***  ROUTINE cice_sbc_in  ***
278      !! ** Purpose: Set coupling fields and pass to CICE
279      !!---------------------------------------------------------------------
280      INTEGER, INTENT(in   ) ::   kt   ! ocean time step
281      INTEGER, INTENT(in   ) ::   ksbc ! surface forcing type
282      !
283      INTEGER  ::   ji, jj, jl                   ! dummy loop indices     
284      REAL(wp), DIMENSION(jpi,jpj) :: ztmp, zpice
285      REAL(wp), DIMENSION(jpi,jpj,ncat) :: ztmpn
286      REAL(wp) ::   zintb, zintn  ! dummy argument
287      !!---------------------------------------------------------------------
288      !
289      IF( kt == nit000 )  THEN
290         IF(lwp) WRITE(numout,*)'cice_sbc_in'
291      ENDIF
292
293      ztmp(:,:)=0.0
294
295! Aggregate ice concentration already set in cice_sbc_out (or cice_sbc_init on
296! the first time-step)
297
298! forced and coupled case
299
300      IF ( (ksbc == jp_flx).OR.(ksbc == jp_purecpl) ) THEN
301
302         ztmpn(:,:,:)=0.0
303
304! x comp of wind stress (CI_1)
305! U point to F point
306         DO jj=1,jpjm1
307            DO ji=1,jpi
308               ztmp(ji,jj) = 0.5 * (  fr_iu(ji,jj) * utau(ji,jj)      &
309                                    + fr_iu(ji,jj+1) * utau(ji,jj+1) ) * fmask(ji,jj,1)
310            ENDDO
311         ENDDO
312         CALL nemo2cice(ztmp,strax,'F', -1. )
313
314! y comp of wind stress (CI_2)
315! V point to F point
316         DO jj=1,jpj
317            DO ji=1,jpim1
318               ztmp(ji,jj) = 0.5 * (  fr_iv(ji,jj) * vtau(ji,jj)      &
319                                    + fr_iv(ji+1,jj) * vtau(ji+1,jj) ) * fmask(ji,jj,1)
320            ENDDO
321         ENDDO
322         CALL nemo2cice(ztmp,stray,'F', -1. )
323
324! Surface downward latent heat flux (CI_5)
325         IF (ksbc == jp_flx) THEN
326            DO jl=1,ncat
327               ztmpn(:,:,jl)=qla_ice(:,:,1)*a_i(:,:,jl)
328            ENDDO
329         ELSE
330! emp_ice is set in sbc_cpl_ice_flx as sublimation-snow
331            qla_ice(:,:,1)= - ( emp_ice(:,:)+sprecip(:,:) ) * rLsub
332! End of temporary code
333            DO jj=1,jpj
334               DO ji=1,jpi
335                  IF (fr_i(ji,jj).eq.0.0) THEN
336                     DO jl=1,ncat
337                        ztmpn(ji,jj,jl)=0.0
338                     ENDDO
339                     ! This will then be conserved in CICE
340                     ztmpn(ji,jj,1)=qla_ice(ji,jj,1)
341                  ELSE
342                     DO jl=1,ncat
343                        ztmpn(ji,jj,jl)=qla_ice(ji,jj,1)*a_i(ji,jj,jl)/fr_i(ji,jj)
344                     ENDDO
345                  ENDIF
346               ENDDO
347            ENDDO
348         ENDIF
349         DO jl=1,ncat
350            CALL nemo2cice(ztmpn(:,:,jl),flatn_f(:,:,jl,:),'T', 1. )
351
352! GBM conductive flux through ice (CI_6)
353!  Convert to GBM
354            IF (ksbc == jp_flx) THEN
355               ztmp(:,:) = botmelt(:,:,jl)*a_i(:,:,jl)
356            ELSE
357               ztmp(:,:) = botmelt(:,:,jl)
358            ENDIF
359            CALL nemo2cice(ztmp,fcondtopn_f(:,:,jl,:),'T', 1. )
360
361! GBM surface heat flux (CI_7)
362!  Convert to GBM
363            IF (ksbc == jp_flx) THEN
364               ztmp(:,:) = (topmelt(:,:,jl)+botmelt(:,:,jl))*a_i(:,:,jl) 
365            ELSE
366               ztmp(:,:) = (topmelt(:,:,jl)+botmelt(:,:,jl))
367            ENDIF
368            CALL nemo2cice(ztmp,fsurfn_f(:,:,jl,:),'T', 1. )
369         ENDDO
370
371      ELSE IF (ksbc == jp_blk) THEN
372
373! Pass bulk forcing fields to CICE (which will calculate heat fluxes etc itself)
374! x comp and y comp of atmosphere surface wind (CICE expects on T points)
375         ztmp(:,:) = wndi_ice(:,:)
376         CALL nemo2cice(ztmp,uatm,'T', -1. )
377         ztmp(:,:) = wndj_ice(:,:)
378         CALL nemo2cice(ztmp,vatm,'T', -1. )
379         ztmp(:,:) = SQRT ( wndi_ice(:,:)**2 + wndj_ice(:,:)**2 )
380         CALL nemo2cice(ztmp,wind,'T', 1. )    ! Wind speed (m/s)
381         ztmp(:,:) = qsr_ice(:,:,1)
382         CALL nemo2cice(ztmp,fsw,'T', 1. )     ! Incoming short-wave (W/m^2)
383         ztmp(:,:) = qlw_ice(:,:,1)
384         CALL nemo2cice(ztmp,flw,'T', 1. )     ! Incoming long-wave (W/m^2)
385         ztmp(:,:) = tatm_ice(:,:)
386         CALL nemo2cice(ztmp,Tair,'T', 1. )    ! Air temperature (K)
387         CALL nemo2cice(ztmp,potT,'T', 1. )    ! Potential temp (K)
388! Following line uses MAX(....) to avoid problems if tatm_ice has unset halo rows 
389         ztmp(:,:) = 101000. / ( 287.04 * MAX(1.0,tatm_ice(:,:)) )   
390                                               ! Constant (101000.) atm pressure assumed
391         CALL nemo2cice(ztmp,rhoa,'T', 1. )    ! Air density (kg/m^3)
392         ztmp(:,:) = qatm_ice(:,:)
393         CALL nemo2cice(ztmp,Qa,'T', 1. )      ! Specific humidity (kg/kg)
394         ztmp(:,:)=10.0
395         CALL nemo2cice(ztmp,zlvl,'T', 1. )    ! Atmos level height (m)
396
397! May want to check all values are physically realistic (as in CICE routine
398! prepare_forcing)?
399
400! Divide shortwave into spectral bands (as in prepare_forcing)
401         ztmp(:,:)=qsr_ice(:,:,1)*frcvdr       ! visible direct
402         CALL nemo2cice(ztmp,swvdr,'T', 1. )             
403         ztmp(:,:)=qsr_ice(:,:,1)*frcvdf       ! visible diffuse
404         CALL nemo2cice(ztmp,swvdf,'T', 1. )             
405         ztmp(:,:)=qsr_ice(:,:,1)*frcidr       ! near IR direct
406         CALL nemo2cice(ztmp,swidr,'T', 1. )
407         ztmp(:,:)=qsr_ice(:,:,1)*frcidf       ! near IR diffuse
408         CALL nemo2cice(ztmp,swidf,'T', 1. )
409
410      ENDIF
411
412! Snowfall
413! Ensure fsnow is positive (as in CICE routine prepare_forcing)
414      IF( iom_use('snowpre') )   CALL iom_put('snowpre',MAX( (1.0-fr_i(:,:))*sprecip(:,:) ,0.0)) !!Joakim edit 
415      ztmp(:,:)=MAX(fr_i(:,:)*sprecip(:,:),0.0) 
416      CALL nemo2cice(ztmp,fsnow,'T', 1. ) 
417
418! Rainfall
419      IF( iom_use('precip') )   CALL iom_put('precip', (1.0-fr_i(:,:))*(tprecip(:,:)-sprecip(:,:)) ) !!Joakim edit
420      ztmp(:,:)=fr_i(:,:)*(tprecip(:,:)-sprecip(:,:))
421      CALL nemo2cice(ztmp,frain,'T', 1. ) 
422
423! Freezing/melting potential
424! Calculated over NEMO leapfrog timestep (hence 2*dt)
425      nfrzmlt(:,:) = rau0 * rcp * e3t_m(:,:) * ( Tocnfrz-sst_m(:,:) ) / ( 2.0*dt )
426
427      ztmp(:,:) = nfrzmlt(:,:)
428      CALL nemo2cice(ztmp,frzmlt,'T', 1. )
429
430! SST  and SSS
431
432      CALL nemo2cice(sst_m,sst,'T', 1. )
433      CALL nemo2cice(sss_m,sss,'T', 1. )
434
435! x comp and y comp of surface ocean current
436! U point to F point
437      DO jj=1,jpjm1
438         DO ji=1,jpi
439            ztmp(ji,jj)=0.5*(ssu_m(ji,jj)+ssu_m(ji,jj+1))*fmask(ji,jj,1)
440         ENDDO
441      ENDDO
442      CALL nemo2cice(ztmp,uocn,'F', -1. )
443
444! V point to F point
445      DO jj=1,jpj
446         DO ji=1,jpim1
447            ztmp(ji,jj)=0.5*(ssv_m(ji,jj)+ssv_m(ji+1,jj))*fmask(ji,jj,1)
448         ENDDO
449      ENDDO
450      CALL nemo2cice(ztmp,vocn,'F', -1. )
451
452      IF( ln_ice_embd ) THEN             !== embedded sea ice: compute representative ice top surface ==!
453          !
454          ! average interpolation coeff as used in dynspg = (1/nn_fsbc) * {SUM[n/nn_fsbc], n=0,nn_fsbc-1}
455          !                                               = (1/nn_fsbc)^2 * {SUM[n], n=0,nn_fsbc-1}
456         zintn = REAL( nn_fsbc - 1 ) / REAL( nn_fsbc ) * 0.5_wp
457          !
458          ! average interpolation coeff as used in dynspg = (1/nn_fsbc) * {SUM[1-n/nn_fsbc], n=0,nn_fsbc-1}
459          !                                               = (1/nn_fsbc)^2 * (nn_fsbc^2 - {SUM[n], n=0,nn_fsbc-1})
460         zintb = REAL( nn_fsbc + 1 ) / REAL( nn_fsbc ) * 0.5_wp
461          !
462         zpice(:,:) = ssh_m(:,:) + (  zintn * snwice_mass(:,:) +  zintb * snwice_mass_b(:,:)  ) * r1_rau0
463          !
464         !
465      ELSE                                    !== non-embedded sea ice: use ocean surface for slope calculation ==!
466         zpice(:,:) = ssh_m(:,:)
467      ENDIF
468
469! x comp and y comp of sea surface slope (on F points)
470! T point to F point
471      DO jj = 1, jpjm1
472         DO ji = 1, jpim1
473            ztmp(ji,jj)=0.5 * (  (zpice(ji+1,jj  )-zpice(ji,jj  )) * r1_e1u(ji,jj  )    &
474               &               + (zpice(ji+1,jj+1)-zpice(ji,jj+1)) * r1_e1u(ji,jj+1)  ) * fmask(ji,jj,1)
475         END DO
476      END DO
477      CALL nemo2cice( ztmp,ss_tltx,'F', -1. )
478
479! T point to F point
480      DO jj = 1, jpjm1
481         DO ji = 1, jpim1
482            ztmp(ji,jj)=0.5 * (  (zpice(ji  ,jj+1)-zpice(ji  ,jj)) * r1_e2v(ji  ,jj)    &
483               &               + (zpice(ji+1,jj+1)-zpice(ji+1,jj)) * r1_e2v(ji+1,jj)  ) *  fmask(ji,jj,1)
484         END DO
485      END DO
486      CALL nemo2cice(ztmp,ss_tlty,'F', -1. )
487      !
488   END SUBROUTINE cice_sbc_in
489
490
491   SUBROUTINE cice_sbc_out( kt, ksbc )
492      !!---------------------------------------------------------------------
493      !!                    ***  ROUTINE cice_sbc_out  ***
494      !! ** Purpose: Get fields from CICE and set surface fields for NEMO
495      !!---------------------------------------------------------------------
496      INTEGER, INTENT( in  ) ::   kt   ! ocean time step
497      INTEGER, INTENT( in  ) ::   ksbc ! surface forcing type
498     
499      INTEGER  ::   ji, jj, jl                 ! dummy loop indices
500      REAL(wp), DIMENSION(jpi,jpj) :: ztmp1, ztmp2
501      !!---------------------------------------------------------------------
502      !
503      IF( kt == nit000 )  THEN
504         IF(lwp) WRITE(numout,*)'cice_sbc_out'
505      ENDIF
506     
507! x comp of ocean-ice stress
508      CALL cice2nemo(strocnx,ztmp1,'F', -1. )
509      ss_iou(:,:)=0.0
510! F point to U point
511      DO jj=2,jpjm1
512         DO ji=2,jpim1
513            ss_iou(ji,jj) = 0.5 * ( ztmp1(ji,jj-1) + ztmp1(ji,jj) ) * umask(ji,jj,1)
514         ENDDO
515      ENDDO
516      CALL lbc_lnk( 'sbcice_cice', ss_iou , 'U', -1. )
517
518! y comp of ocean-ice stress
519      CALL cice2nemo(strocny,ztmp1,'F', -1. )
520      ss_iov(:,:)=0.0
521! F point to V point
522
523      DO jj=1,jpjm1
524         DO ji=2,jpim1
525            ss_iov(ji,jj) = 0.5 * ( ztmp1(ji-1,jj) + ztmp1(ji,jj) ) * vmask(ji,jj,1)
526         ENDDO
527      ENDDO
528      CALL lbc_lnk( 'sbcice_cice', ss_iov , 'V', -1. )
529
530! x and y comps of surface stress
531! Combine wind stress and ocean-ice stress
532! [Note that fr_iu hasn't yet been updated, so still from start of CICE timestep]
533! strocnx and strocny already weighted by ice fraction in CICE so not done here
534
535      utau(:,:)=(1.0-fr_iu(:,:))*utau(:,:)-ss_iou(:,:)
536      vtau(:,:)=(1.0-fr_iv(:,:))*vtau(:,:)-ss_iov(:,:)     
537 
538! Also need ice/ocean stress on T points so that taum can be updated
539! This interpolation is already done in CICE so best to use those values
540      CALL cice2nemo(strocnxT,ztmp1,'T',-1.) 
541      CALL cice2nemo(strocnyT,ztmp2,'T',-1.) 
542 
543! Update taum with modulus of ice-ocean stress
544! strocnxT and strocnyT are not weighted by ice fraction in CICE so must be done here
545taum(:,:)=(1.0-fr_i(:,:))*taum(:,:)+fr_i(:,:)*SQRT(ztmp1*ztmp1 + ztmp2*ztmp2) 
546
547! Freshwater fluxes
548
549      IF (ksbc == jp_flx) THEN
550! Note that emp from the forcing files is evap*(1-aice)-(tprecip-aice*sprecip)
551! What we want here is evap*(1-aice)-tprecip*(1-aice) hence manipulation below
552! Not ideal since aice won't be the same as in the atmosphere. 
553! Better to use evap and tprecip? (but for now don't read in evap in this case)
554         emp(:,:)  = emp(:,:)+fr_i(:,:)*(tprecip(:,:)-sprecip(:,:))
555      ELSE IF (ksbc == jp_blk) THEN
556         emp(:,:)  = (1.0-fr_i(:,:))*emp(:,:)       
557      ELSE IF (ksbc == jp_purecpl) THEN
558! emp_tot is set in sbc_cpl_ice_flx (called from cice_sbc_in above)
559! This is currently as required with the coupling fields from the UM atmosphere
560         emp(:,:) = emp_tot(:,:)+tprecip(:,:)*fr_i(:,:) 
561      ENDIF
562
563#if defined key_cice4
564      CALL cice2nemo(fresh_gbm,ztmp1,'T', 1. )
565      CALL cice2nemo(fsalt_gbm,ztmp2,'T', 1. )
566#else
567      CALL cice2nemo(fresh_ai,ztmp1,'T', 1. )
568      CALL cice2nemo(fsalt_ai,ztmp2,'T', 1. )
569#endif
570
571! Check to avoid unphysical expression when ice is forming (ztmp1 negative)
572! Otherwise we are effectively allowing ice of higher salinity than the ocean to form
573! which has to be compensated for by the ocean salinity potentially going negative
574! This check breaks conservation but seems reasonable until we have prognostic ice salinity
575! Note the 1000.0 below is to convert from kg salt to g salt (needed for PSU)
576      WHERE (ztmp1(:,:).lt.0.0) ztmp2(:,:)=MAX(ztmp2(:,:),ztmp1(:,:)*sss_m(:,:)/1000.0)
577      sfx(:,:)=ztmp2(:,:)*1000.0
578      emp(:,:)=emp(:,:)-ztmp1(:,:)
579      fmmflx(:,:) = ztmp1(:,:) !!Joakim edit
580     
581      CALL lbc_lnk_multi( 'sbcice_cice', emp , 'T', 1., sfx , 'T', 1. )
582
583! Solar penetrative radiation and non solar surface heat flux
584
585! Scale qsr and qns according to ice fraction (bulk formulae only)
586
587      IF (ksbc == jp_blk) THEN
588         qsr(:,:)=qsr(:,:)*(1.0-fr_i(:,:))
589         qns(:,:)=qns(:,:)*(1.0-fr_i(:,:))
590      ENDIF
591! Take into account snow melting except for fully coupled when already in qns_tot
592      IF (ksbc == jp_purecpl) THEN
593         qsr(:,:)= qsr_tot(:,:)
594         qns(:,:)= qns_tot(:,:)
595      ELSE
596         qns(:,:)= qns(:,:)-sprecip(:,:)*Lfresh*(1.0-fr_i(:,:))
597      ENDIF
598
599! Now add in ice / snow related terms
600! [fswthru will be zero unless running with calc_Tsfc=T in CICE]
601#if defined key_cice4
602      CALL cice2nemo(fswthru_gbm,ztmp1,'T', 1. )
603#else
604      CALL cice2nemo(fswthru_ai,ztmp1,'T', 1. )
605#endif
606      qsr(:,:)=qsr(:,:)+ztmp1(:,:)
607      CALL lbc_lnk( 'sbcice_cice', qsr , 'T', 1. )
608
609      DO jj=1,jpj
610         DO ji=1,jpi
611            nfrzmlt(ji,jj)=MAX(nfrzmlt(ji,jj),0.0)
612         ENDDO
613      ENDDO
614
615#if defined key_cice4
616      CALL cice2nemo(fhocn_gbm,ztmp1,'T', 1. )
617#else
618      CALL cice2nemo(fhocn_ai,ztmp1,'T', 1. )
619#endif
620      qns(:,:)=qns(:,:)+nfrzmlt(:,:)+ztmp1(:,:)
621
622      CALL lbc_lnk( 'sbcice_cice', qns , 'T', 1. )
623
624! Prepare for the following CICE time-step
625
626      CALL cice2nemo(aice,fr_i,'T', 1. )
627      IF ( (ksbc == jp_flx).OR.(ksbc == jp_purecpl) ) THEN
628         DO jl=1,ncat
629            CALL cice2nemo(aicen(:,:,jl,:),a_i(:,:,jl), 'T', 1. )
630         ENDDO
631      ENDIF
632
633! T point to U point
634! T point to V point
635      DO jj=1,jpjm1
636         DO ji=1,jpim1
637            fr_iu(ji,jj)=0.5*(fr_i(ji,jj)+fr_i(ji+1,jj))*umask(ji,jj,1)
638            fr_iv(ji,jj)=0.5*(fr_i(ji,jj)+fr_i(ji,jj+1))*vmask(ji,jj,1)
639         ENDDO
640      ENDDO
641
642      CALL lbc_lnk_multi( 'sbcice_cice', fr_iu , 'U', 1., fr_iv , 'V', 1. )
643
644      ! set the snow+ice mass
645      CALL cice2nemo(vsno(:,:,:),ztmp1,'T', 1. )
646      CALL cice2nemo(vice(:,:,:),ztmp2,'T', 1. )
647      snwice_mass  (:,:) = ( rhos * ztmp1(:,:) + rhoi * ztmp2(:,:)  )
648      snwice_mass_b(:,:) = snwice_mass(:,:)
649      snwice_fmass (:,:) = ( snwice_mass(:,:) - snwice_mass_b(:,:) ) / dt
650      !
651   END SUBROUTINE cice_sbc_out
652
653
654   SUBROUTINE cice_sbc_hadgam( kt )
655      !!---------------------------------------------------------------------
656      !!                    ***  ROUTINE cice_sbc_hadgam  ***
657      !! ** Purpose: Prepare fields needed to pass to HadGAM3 atmosphere
658      !!
659      !!
660      !!---------------------------------------------------------------------
661      INTEGER, INTENT( in  ) ::   kt   ! ocean time step
662      !!
663      INTEGER  ::   jl                        ! dummy loop index
664      INTEGER  ::   ierror
665      !!---------------------------------------------------------------------
666      !
667      IF( kt == nit000 )  THEN
668         IF(lwp) WRITE(numout,*)'cice_sbc_hadgam'
669         IF( sbc_cpl_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'sbc_cpl_alloc : unable to allocate arrays' )
670      ENDIF
671
672      !                                         ! =========================== !
673      !                                         !   Prepare Coupling fields   !
674      !                                         ! =========================== !
675      !
676      ! x and y comp of ice velocity
677      !
678      CALL cice2nemo(uvel,u_ice,'F', -1. )
679      CALL cice2nemo(vvel,v_ice,'F', -1. )
680      !
681      ! Ice concentration (CO_1) = a_i calculated at end of cice_sbc_out 
682      !
683      ! Snow and ice thicknesses (CO_2 and CO_3)
684      !
685      DO jl = 1, ncat
686         CALL cice2nemo( vsnon(:,:,jl,:), h_s(:,:,jl),'T', 1. )
687         CALL cice2nemo( vicen(:,:,jl,:), h_i(:,:,jl),'T', 1. )
688      END DO
689      !
690   END SUBROUTINE cice_sbc_hadgam
691
692
693   SUBROUTINE cice_sbc_final
694      !!---------------------------------------------------------------------
695      !!                    ***  ROUTINE cice_sbc_final  ***
696      !! ** Purpose: Finalize CICE
697      !!---------------------------------------------------------------------
698      !
699      IF(lwp) WRITE(numout,*)'cice_sbc_final'
700      !
701      CALL CICE_Finalize
702      !
703   END SUBROUTINE cice_sbc_final
704
705
706   SUBROUTINE cice_sbc_force (kt)
707      !!---------------------------------------------------------------------
708      !!                    ***  ROUTINE cice_sbc_force  ***
709      !! ** Purpose : Provide CICE forcing from files
710      !!
711      !!---------------------------------------------------------------------
712      !! ** Method  :   READ monthly flux file in NetCDF files
713      !!     
714      !!  snowfall   
715      !!  rainfall   
716      !!  sublimation rate   
717      !!  topmelt (category)
718      !!  botmelt (category)
719      !!
720      !! History :
721      !!----------------------------------------------------------------------
722      USE iom
723      !!
724      INTEGER, INTENT( in  ) ::   kt   ! ocean time step
725      !!
726      INTEGER  ::   ierror             ! return error code
727      INTEGER  ::   ifpr               ! dummy loop index
728      !!
729      CHARACTER(len=100) ::  cn_dir                            !   Root directory for location of CICE forcing files
730      TYPE(FLD_N), DIMENSION(jpfld) ::   slf_i                 ! array of namelist informations on the fields to read
731      TYPE(FLD_N) ::   sn_snow, sn_rain, sn_sblm               ! informations about the fields to be read
732      TYPE(FLD_N) ::   sn_top1, sn_top2, sn_top3, sn_top4, sn_top5
733      TYPE(FLD_N) ::   sn_bot1, sn_bot2, sn_bot3, sn_bot4, sn_bot5 
734      !!
735      NAMELIST/namsbc_cice/ cn_dir, sn_snow, sn_rain, sn_sblm,   &
736         &                          sn_top1, sn_top2, sn_top3, sn_top4, sn_top5,  &
737         &                          sn_bot1, sn_bot2, sn_bot3, sn_bot4, sn_bot5
738      INTEGER :: ios
739      !!---------------------------------------------------------------------
740
741      !                                         ! ====================== !
742      IF( kt == nit000 ) THEN                   !  First call kt=nit000  !
743         !                                      ! ====================== !
744         ! namsbc_cice is not yet in the reference namelist
745         ! set file information (default values)
746         cn_dir = './'       ! directory in which the model is executed
747
748         ! (NB: frequency positive => hours, negative => months)
749         !            !    file          ! frequency !  variable    ! time intep !  clim   ! 'yearly' or ! weights  ! rotation   ! landmask
750         !            !    name          !  (hours)  !   name       !   (T/F)    !  (T/F)  !  'monthly'  ! filename ! pairs      ! file
751         sn_snow = FLD_N( 'snowfall_1m'  ,    -1.    ,  'snowfall'  ,  .true.    , .true.  ,  ' yearly'  , ''       , ''         ,  ''    ) 
752         sn_rain = FLD_N( 'rainfall_1m'  ,    -1.    ,  'rainfall'  ,  .true.    , .true.  ,  ' yearly'  , ''       , ''         ,  ''    ) 
753         sn_sblm = FLD_N( 'sublim_1m'    ,    -1.    ,  'sublim'    ,  .true.    , .true.  ,  ' yearly'  , ''       , ''         ,  ''    )
754         sn_top1 = FLD_N( 'topmeltn1_1m' ,    -1.    ,  'topmeltn1' ,  .true.    , .true.  ,  ' yearly'  , ''       , ''         ,  ''    )
755         sn_top2 = FLD_N( 'topmeltn2_1m' ,    -1.    ,  'topmeltn2' ,  .true.    , .true.  ,  ' yearly'  , ''       , ''         ,  ''    )
756         sn_top3 = FLD_N( 'topmeltn3_1m' ,    -1.    ,  'topmeltn3' ,  .true.    , .true.  ,  ' yearly'  , ''       , ''         ,  ''    )
757         sn_top4 = FLD_N( 'topmeltn4_1m' ,    -1.    ,  'topmeltn4' ,  .true.    , .true.  ,  ' yearly'  , ''       , ''         ,  ''    )
758         sn_top5 = FLD_N( 'topmeltn5_1m' ,    -1.    ,  'topmeltn5' ,  .true.    , .true.  ,  ' yearly'  , ''       , ''         ,  ''    )
759         sn_bot1 = FLD_N( 'botmeltn1_1m' ,    -1.    ,  'botmeltn1' ,  .true.    , .true.  ,  ' yearly'  , ''       , ''         ,  ''    )
760         sn_bot2 = FLD_N( 'botmeltn2_1m' ,    -1.    ,  'botmeltn2' ,  .true.    , .true.  ,  ' yearly'  , ''       , ''         ,  ''    )
761         sn_bot3 = FLD_N( 'botmeltn3_1m' ,    -1.    ,  'botmeltn3' ,  .true.    , .true.  ,  ' yearly'  , ''       , ''         ,  ''    )
762         sn_bot4 = FLD_N( 'botmeltn4_1m' ,    -1.    ,  'botmeltn4' ,  .true.    , .true.  ,  ' yearly'  , ''       , ''         ,  ''    )
763         sn_bot5 = FLD_N( 'botmeltn5_1m' ,    -1.    ,  'botmeltn5' ,  .true.    , .true.  ,  ' yearly'  , ''       , ''         ,  ''    )
764
765         READ  ( numnam_ref, namsbc_cice, IOSTAT = ios, ERR = 901)
766901      IF( ios /= 0 )   CALL ctl_nam ( ios , 'namsbc_cice in reference namelist' )
767
768         READ  ( numnam_cfg, namsbc_cice, IOSTAT = ios, ERR = 902 )
769902      IF( ios >  0 )   CALL ctl_nam ( ios , 'namsbc_cice in configuration namelist' )
770         IF(lwm) WRITE ( numond, namsbc_cice )
771
772         ! store namelist information in an array
773         slf_i(jp_snow) = sn_snow   ;   slf_i(jp_rain) = sn_rain   ;   slf_i(jp_sblm) = sn_sblm
774         slf_i(jp_top1) = sn_top1   ;   slf_i(jp_top2) = sn_top2   ;   slf_i(jp_top3) = sn_top3
775         slf_i(jp_top4) = sn_top4   ;   slf_i(jp_top5) = sn_top5   ;   slf_i(jp_bot1) = sn_bot1
776         slf_i(jp_bot2) = sn_bot2   ;   slf_i(jp_bot3) = sn_bot3   ;   slf_i(jp_bot4) = sn_bot4
777         slf_i(jp_bot5) = sn_bot5
778         
779         ! set sf structure
780         ALLOCATE( sf(jpfld), STAT=ierror )
781         IF( ierror > 0 ) THEN
782            CALL ctl_stop( 'cice_sbc_force: unable to allocate sf structure' )   ;   RETURN
783         ENDIF
784
785         DO ifpr= 1, jpfld
786            ALLOCATE( sf(ifpr)%fnow(jpi,jpj,1) )
787            ALLOCATE( sf(ifpr)%fdta(jpi,jpj,1,2) )
788         END DO
789
790         ! fill sf with slf_i and control print
791         CALL fld_fill( sf, slf_i, cn_dir, 'cice_sbc_force', 'flux formulation for CICE', 'namsbc_cice' )
792         !
793      ENDIF
794
795      CALL fld_read( kt, nn_fsbc, sf )           ! Read input fields and provides the
796      !                                          ! input fields at the current time-step
797
798      ! set the fluxes from read fields
799      sprecip(:,:) = sf(jp_snow)%fnow(:,:,1)
800      tprecip(:,:) = sf(jp_snow)%fnow(:,:,1)+sf(jp_rain)%fnow(:,:,1)
801! May be better to do this conversion somewhere else
802      qla_ice(:,:,1) = -rLsub*sf(jp_sblm)%fnow(:,:,1)
803      topmelt(:,:,1) = sf(jp_top1)%fnow(:,:,1)
804      topmelt(:,:,2) = sf(jp_top2)%fnow(:,:,1)
805      topmelt(:,:,3) = sf(jp_top3)%fnow(:,:,1)
806      topmelt(:,:,4) = sf(jp_top4)%fnow(:,:,1)
807      topmelt(:,:,5) = sf(jp_top5)%fnow(:,:,1)
808      botmelt(:,:,1) = sf(jp_bot1)%fnow(:,:,1)
809      botmelt(:,:,2) = sf(jp_bot2)%fnow(:,:,1)
810      botmelt(:,:,3) = sf(jp_bot3)%fnow(:,:,1)
811      botmelt(:,:,4) = sf(jp_bot4)%fnow(:,:,1)
812      botmelt(:,:,5) = sf(jp_bot5)%fnow(:,:,1)
813
814      ! control print (if less than 100 time-step asked)
815      IF( nitend-nit000 <= 100 .AND. lwp ) THEN
816         WRITE(numout,*) 
817         WRITE(numout,*) '        read forcing fluxes for CICE OK'
818         CALL FLUSH(numout)
819      ENDIF
820
821   END SUBROUTINE cice_sbc_force
822
823   SUBROUTINE nemo2cice( pn, pc, cd_type, psgn)
824      !!---------------------------------------------------------------------
825      !!                    ***  ROUTINE nemo2cice  ***
826      !! ** Purpose :   Transfer field in NEMO array to field in CICE array. 
827#if defined key_nemocice_decomp
828      !!             
829      !!                NEMO and CICE PE sub domains are identical, hence
830      !!                there is no need to gather or scatter data from
831      !!                one PE configuration to another.
832#else
833      !!                Automatically gather/scatter between
834      !!                different processors and blocks
835      !! ** Method :    A. Ensure all haloes are filled in NEMO field (pn)
836      !!                B. Gather pn into global array (png)
837      !!                C. Map png into CICE global array (pcg)
838      !!                D. Scatter pcg to CICE blocks (pc) + update haloes 
839#endif
840      !!---------------------------------------------------------------------
841      CHARACTER(len=1), INTENT( in ) ::   &
842          cd_type       ! nature of pn grid-point
843          !             !   = T or F gridpoints
844      REAL(wp), INTENT( in ) ::   &
845          psgn          ! control of the sign change
846          !             !   =-1 , the sign is modified following the type of b.c. used
847          !             !   = 1 , no sign change
848      REAL(wp), DIMENSION(jpi,jpj) :: pn
849#if !defined key_nemocice_decomp
850      REAL(wp), DIMENSION(jpiglo,jpjglo) :: png2
851      REAL (kind=dbl_kind), dimension(nx_global,ny_global) :: pcg
852#endif
853      REAL (kind=dbl_kind), dimension(nx_block,ny_block,max_blocks) :: pc
854      INTEGER (int_kind) :: &
855         field_type,        &! id for type of field (scalar, vector, angle)
856         grid_loc            ! id for location on horizontal grid
857                            !  (center, NEcorner, Nface, Eface)
858
859      INTEGER  ::   ji, jj, jn                      ! dummy loop indices
860      !!---------------------------------------------------------------------
861
862!     A. Ensure all haloes are filled in NEMO field (pn)
863
864      CALL lbc_lnk( 'sbcice_cice', pn , cd_type, psgn )
865
866#if defined key_nemocice_decomp
867
868      ! Copy local domain data from NEMO to CICE field
869      pc(:,:,1)=0.0
870      DO jj=2,ny_block-1
871         DO ji=2,nx_block-1
872            pc(ji,jj,1)=pn(ji-1+ji_off,jj-1+jj_off)
873         ENDDO
874      ENDDO
875
876#else
877
878!     B. Gather pn into global array (png)
879
880      IF ( jpnij > 1) THEN
881         CALL mppsync
882         CALL mppgather (pn,0,png) 
883         CALL mppsync
884      ELSE
885         png(:,:,1)=pn(:,:)
886      ENDIF
887
888!     C. Map png into CICE global array (pcg)
889
890! Need to make sure this is robust to changes in NEMO halo rows....
891! (may be OK but not 100% sure)
892
893      IF (nproc==0) THEN     
894!        pcg(:,:)=0.0
895         DO jn=1,jpnij
896            DO jj=nldjt(jn),nlejt(jn)
897               DO ji=nldit(jn),nleit(jn)
898                  png2(ji+nimppt(jn)-1,jj+njmppt(jn)-1)=png(ji,jj,jn)
899               ENDDO
900            ENDDO
901         ENDDO
902         DO jj=1,ny_global
903            DO ji=1,nx_global
904               pcg(ji,jj)=png2(ji+ji_off,jj+jj_off)
905            ENDDO
906         ENDDO
907      ENDIF
908
909#endif
910
911      SELECT CASE ( cd_type )
912         CASE ( 'T' )
913            grid_loc=field_loc_center
914         CASE ( 'F' )                             
915            grid_loc=field_loc_NEcorner
916      END SELECT
917
918      SELECT CASE ( NINT(psgn) )
919         CASE ( -1 )
920            field_type=field_type_vector
921         CASE ( 1 )                             
922            field_type=field_type_scalar
923      END SELECT
924
925#if defined key_nemocice_decomp
926      ! Ensure CICE halos are up to date
927      CALL ice_HaloUpdate (pc, halo_info, grid_loc, field_type)
928#else
929!     D. Scatter pcg to CICE blocks (pc) + update halos
930      CALL scatter_global(pc, pcg, 0, distrb_info, grid_loc, field_type)
931#endif
932
933   END SUBROUTINE nemo2cice
934
935   SUBROUTINE cice2nemo ( pc, pn, cd_type, psgn )
936      !!---------------------------------------------------------------------
937      !!                    ***  ROUTINE cice2nemo  ***
938      !! ** Purpose :   Transfer field in CICE array to field in NEMO array.
939#if defined key_nemocice_decomp
940      !!             
941      !!                NEMO and CICE PE sub domains are identical, hence
942      !!                there is no need to gather or scatter data from
943      !!                one PE configuration to another.
944#else 
945      !!                Automatically deal with scatter/gather between
946      !!                different processors and blocks
947      !! ** Method :    A. Gather CICE blocks (pc) into global array (pcg)
948      !!                B. Map pcg into NEMO global array (png)
949      !!                C. Scatter png into NEMO field (pn) for each processor
950      !!                D. Ensure all haloes are filled in pn
951#endif
952      !!---------------------------------------------------------------------
953
954      CHARACTER(len=1), INTENT( in ) ::   &
955          cd_type       ! nature of pn grid-point
956          !             !   = T or F gridpoints
957      REAL(wp), INTENT( in ) ::   &
958          psgn          ! control of the sign change
959          !             !   =-1 , the sign is modified following the type of b.c. used
960          !             !   = 1 , no sign change
961      REAL(wp), DIMENSION(jpi,jpj) :: pn
962
963#if defined key_nemocice_decomp
964      INTEGER (int_kind) :: &
965         field_type,        & ! id for type of field (scalar, vector, angle)
966         grid_loc             ! id for location on horizontal grid
967                              ! (center, NEcorner, Nface, Eface)
968#else
969      REAL (kind=dbl_kind), dimension(nx_global,ny_global) :: pcg
970#endif
971
972      REAL (kind=dbl_kind), dimension(nx_block,ny_block,max_blocks) :: pc
973
974      INTEGER  ::   ji, jj, jn                      ! dummy loop indices
975
976
977#if defined key_nemocice_decomp
978
979      SELECT CASE ( cd_type )
980         CASE ( 'T' )
981            grid_loc=field_loc_center
982         CASE ( 'F' )                             
983            grid_loc=field_loc_NEcorner
984      END SELECT
985
986      SELECT CASE ( NINT(psgn) )
987         CASE ( -1 )
988            field_type=field_type_vector
989         CASE ( 1 )                             
990            field_type=field_type_scalar
991      END SELECT
992
993      CALL ice_HaloUpdate (pc, halo_info, grid_loc, field_type)
994
995
996      pn(:,:)=0.0
997      DO jj=1,jpjm1
998         DO ji=1,jpim1
999            pn(ji,jj)=pc(ji+1-ji_off,jj+1-jj_off,1)
1000         ENDDO
1001      ENDDO
1002
1003#else
1004
1005!      A. Gather CICE blocks (pc) into global array (pcg)
1006
1007      CALL gather_global(pcg, pc, 0, distrb_info)
1008
1009!     B. Map pcg into NEMO global array (png)
1010
1011! Need to make sure this is robust to changes in NEMO halo rows....
1012! (may be OK but not spent much time thinking about it)
1013! Note that non-existent pcg elements may be used below, but
1014! the lbclnk call on pn will replace these with sensible values
1015
1016      IF (nproc==0) THEN
1017         png(:,:,:)=0.0
1018         DO jn=1,jpnij
1019            DO jj=nldjt(jn),nlejt(jn)
1020               DO ji=nldit(jn),nleit(jn)
1021                  png(ji,jj,jn)=pcg(ji+nimppt(jn)-1-ji_off,jj+njmppt(jn)-1-jj_off)
1022               ENDDO
1023            ENDDO
1024         ENDDO
1025      ENDIF
1026
1027!     C. Scatter png into NEMO field (pn) for each processor
1028
1029      IF ( jpnij > 1) THEN
1030         CALL mppsync
1031         CALL mppscatter (png,0,pn) 
1032         CALL mppsync
1033      ELSE
1034         pn(:,:)=png(:,:,1)
1035      ENDIF
1036
1037#endif
1038
1039!     D. Ensure all haloes are filled in pn
1040
1041      CALL lbc_lnk( 'sbcice_cice', pn , cd_type, psgn )
1042
1043   END SUBROUTINE cice2nemo
1044
1045#else
1046   !!----------------------------------------------------------------------
1047   !!   Default option           Dummy module         NO CICE sea-ice model
1048   !!----------------------------------------------------------------------
1049CONTAINS
1050
1051   SUBROUTINE sbc_ice_cice ( kt, ksbc )     ! Dummy routine
1052      IMPLICIT NONE
1053      INTEGER, INTENT( in ) :: kt, ksbc
1054      WRITE(*,*) 'sbc_ice_cice: You should not have seen this print! error?', kt
1055   END SUBROUTINE sbc_ice_cice
1056
1057   SUBROUTINE cice_sbc_init (ksbc, Kbb, Kmm)    ! Dummy routine
1058      IMPLICIT NONE
1059      INTEGER, INTENT( in ) :: ksbc
1060      INTEGER, INTENT( in ) :: Kbb, Kmm
1061      WRITE(*,*) 'cice_sbc_init: You should not have seen this print! error?', ksbc
1062   END SUBROUTINE cice_sbc_init
1063
1064   SUBROUTINE cice_sbc_final     ! Dummy routine
1065      IMPLICIT NONE
1066      WRITE(*,*) 'cice_sbc_final: You should not have seen this print! error?'
1067   END SUBROUTINE cice_sbc_final
1068
1069#endif
1070
1071   !!======================================================================
1072END MODULE sbcice_cice
Note: See TracBrowser for help on using the repository browser.