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.
trasbc.F90 in branches/UKMO/dev_isf_flx_UKESM_r9321/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/UKMO/dev_isf_flx_UKESM_r9321/NEMOGCM/NEMO/OPA_SRC/TRA/trasbc.F90 @ 9855

Last change on this file since 9855 was 9855, checked in by mathiot, 6 years ago

update sbcisf from trunk

File size: 16.7 KB
Line 
1MODULE trasbc
2   !!==============================================================================
3   !!                       ***  MODULE  trasbc  ***
4   !! Ocean active tracers:  surface boundary condition
5   !!==============================================================================
6   !! History :  OPA  !  1998-10  (G. Madec, G. Roullet, M. Imbard)  Original code
7   !!            8.2  !  2001-02  (D. Ludicone)  sea ice and free surface
8   !!  NEMO      1.0  !  2002-06  (G. Madec)  F90: Free form and module
9   !!            3.3  !  2010-04  (M. Leclair, G. Madec)  Forcing averaged over 2 time steps
10   !!             -   !  2010-09  (C. Ethe, G. Madec) Merge TRA-TRC
11   !!            3.6  !  2014-11  (P. Mathiot) isf melting forcing
12   !!----------------------------------------------------------------------
13
14   !!----------------------------------------------------------------------
15   !!   tra_sbc      : update the tracer trend at ocean surface
16   !!----------------------------------------------------------------------
17   USE oce             ! ocean dynamics and active tracers
18   USE sbc_oce         ! surface boundary condition: ocean
19   USE dom_oce         ! ocean space domain variables
20   USE phycst          ! physical constant
21   USE sbcmod          ! ln_rnf 
22   USE sbcrnf          ! River runoff 
23   USE sbcisf          ! Ice shelf   
24   USE traqsr          ! solar radiation penetration
25   USE trd_oce         ! trends: ocean variables
26   USE trdtra          ! trends manager: tracers
27   !
28   USE in_out_manager  ! I/O manager
29   USE prtctl          ! Print control
30   USE iom
31   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
32   USE wrk_nemo        ! Memory Allocation
33   USE timing          ! Timing
34   USE eosbn2
35#if defined key_asminc   
36   USE asminc          ! Assimilation increment
37#endif
38
39   IMPLICIT NONE
40   PRIVATE
41
42   PUBLIC   tra_sbc    ! routine called by step.F90
43
44   !! * Substitutions
45#  include "domzgr_substitute.h90"
46#  include "vectopt_loop_substitute.h90"
47   !!----------------------------------------------------------------------
48   !! NEMO/OPA 3.7 , NEMO Consortium (2014)
49   !! $Id$
50   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
51   !!----------------------------------------------------------------------
52CONTAINS
53
54   SUBROUTINE tra_sbc ( kt )
55      !!----------------------------------------------------------------------
56      !!                  ***  ROUTINE tra_sbc  ***
57      !!                   
58      !! ** Purpose :   Compute the tracer surface boundary condition trend of
59      !!      (flux through the interface, concentration/dilution effect)
60      !!      and add it to the general trend of tracer equations.
61      !!
62      !! ** Method :
63      !!      Following Roullet and Madec (2000), the air-sea flux can be divided
64      !!      into three effects: (1) Fext, external forcing;
65      !!      (2) Fwi, concentration/dilution effect due to water exchanged
66      !!         at the surface by evaporation, precipitations and runoff (E-P-R);
67      !!      (3) Fwe, tracer carried with the water that is exchanged.
68      !!            - salinity    : salt flux only due to freezing/melting
69      !!            sa = sa +  sfx / rau0 / e3t  for k=1
70      !!
71      !!      Fext, flux through the air-sea interface for temperature and salt:
72      !!            - temperature : heat flux q (w/m2). If penetrative solar
73      !!         radiation q is only the non solar part of the heat flux, the
74      !!         solar part is added in traqsr.F routine.
75      !!            ta = ta + q /(rau0 rcp e3t)  for k=1
76      !!            - salinity    : no salt flux
77      !!
78      !!      The formulation for Fwb and Fwi vary according to the free
79      !!      surface formulation (linear or variable volume).
80      !!      * Linear free surface
81      !!            The surface freshwater flux modifies the ocean volume
82      !!         and thus the concentration of a tracer and the temperature.
83      !!         First order of the effect of surface freshwater exchange
84      !!         for salinity, it can be neglected on temperature (especially
85      !!         as the temperature of precipitations and runoffs is usually
86      !!         unknown).
87      !!            - temperature : we assume that the temperature of both
88      !!         precipitations and runoffs is equal to the SST, thus there
89      !!         is no additional flux since in this case, the concentration
90      !!         dilution effect is balanced by the net heat flux associated
91      !!         to the freshwater exchange (Fwe+Fwi=0):
92      !!            (Tp P - Te E) + SST (P-E) = 0 when Tp=Te=SST
93      !!            - salinity    : evaporation, precipitation and runoff
94      !!         water has a zero salinity  but there is a salt flux due to
95      !!         freezing/melting, thus:
96      !!            sa = sa + emp * sn / rau0 / e3t   for k=1
97      !!                    + sfx    / rau0 / e3t
98      !!         where emp, the surface freshwater budget (evaporation minus
99      !!         precipitation minus runoff) given in kg/m2/s is divided
100      !!         by rau0 (density of sea water) to obtain m/s.   
101      !!         Note: even though Fwe does not appear explicitly for
102      !!         temperature in this routine, the heat carried by the water
103      !!         exchanged through the surface is part of the total heat flux
104      !!         forcing and must be taken into account in the global heat
105      !!         balance).
106      !!      * nonlinear free surface (variable volume, lk_vvl)
107      !!         contrary to the linear free surface case, Fwi is properly
108      !!         taken into account by using the true layer thicknesses to       
109      !!         calculate tracer content and advection. There is no need to
110      !!         deal with it in this routine.
111      !!           - temperature: Fwe=SST (P-E+R) is added to Fext.
112      !!           - salinity:  Fwe = 0, there is no surface flux of salt.
113      !!
114      !! ** Action  : - Update the 1st level of (ta,sa) with the trend associated
115      !!                with the tracer surface boundary condition
116      !!              - send trends to trdtra module (l_trdtra=T)
117      !!----------------------------------------------------------------------
118      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
119      !!
120      INTEGER  ::   ji, jj, jk, jn           ! dummy loop indices 
121      INTEGER  ::   ikt, ikb 
122      REAL(wp) ::   zfact, z1_e3t, zdep
123      REAL(wp) ::   zt_frz, zpress
124      REAL(wp), POINTER, DIMENSION(:,:,:) ::  ztrdt, ztrds
125      !!----------------------------------------------------------------------
126      !
127      IF( nn_timing == 1 )  CALL timing_start('tra_sbc')
128      !
129      IF( kt == nit000 ) THEN
130         IF(lwp) WRITE(numout,*)
131         IF(lwp) WRITE(numout,*) 'tra_sbc : TRAcer Surface Boundary Condition'
132         IF(lwp) WRITE(numout,*) '~~~~~~~ '
133      ENDIF
134
135      IF( l_trdtra ) THEN                    !* Save ta and sa trends
136         CALL wrk_alloc( jpi, jpj, jpk, ztrdt, ztrds ) 
137         ztrdt(:,:,:) = tsa(:,:,:,jp_tem)
138         ztrds(:,:,:) = tsa(:,:,:,jp_sal)
139      ENDIF
140
141!!gm      IF( .NOT.ln_traqsr )   qsr(:,:) = 0.e0   ! no solar radiation penetration
142      IF( .NOT.ln_traqsr ) THEN     ! no solar radiation penetration
143         qns(:,:) = qns(:,:) + qsr(:,:)      ! total heat flux in qns
144         qsr(:,:) = 0.e0                     ! qsr set to zero
145      ENDIF
146
147      !----------------------------------------
148      !        EMP, SFX and QNS effects
149      !----------------------------------------
150      !                                          Set before sbc tracer content fields
151      !                                          ************************************
152      IF( kt == nit000 ) THEN                      ! Set the forcing field at nit000 - 1
153         !                                         ! -----------------------------------
154         IF( ln_rstart .AND.    &                     ! Restart: read in restart file
155              & iom_varid( numror, 'sbc_hc_b', ldstop = .FALSE. ) > 0 ) THEN
156            IF(lwp) WRITE(numout,*) '          nit000-1 surface tracer content forcing fields red in the restart file'
157            zfact = 0.5_wp
158            CALL iom_get( numror, jpdom_autoglo, 'sbc_hc_b', sbc_tsc_b(:,:,jp_tem) )   ! before heat content sbc trend
159            CALL iom_get( numror, jpdom_autoglo, 'sbc_sc_b', sbc_tsc_b(:,:,jp_sal) )   ! before salt content sbc trend
160         ELSE                                         ! No restart or restart not found: Euler forward time stepping
161            zfact = 1._wp
162            sbc_tsc_b(:,:,:) = 0._wp
163         ENDIF
164      ELSE                                         ! Swap of forcing fields
165         !                                         ! ----------------------
166         zfact = 0.5_wp
167         sbc_tsc_b(:,:,:) = sbc_tsc(:,:,:)
168      ENDIF
169      !                                          Compute now sbc tracer content fields
170      !                                          *************************************
171
172                                                   ! Concentration dilution effect on (t,s) due to 
173                                                   ! evaporation, precipitation and qns, but not river runoff
174                                               
175      IF( lk_vvl ) THEN                            ! Variable Volume case  ==>> heat content of mass flux is in qns
176         DO jj = 1, jpj
177            DO ji = 1, jpi 
178               sbc_tsc(ji,jj,jp_tem) = r1_rau0_rcp * qns(ji,jj)                              ! non solar heat flux
179               sbc_tsc(ji,jj,jp_sal) = r1_rau0     * sfx(ji,jj)                              ! salt flux due to freezing/melting
180            END DO
181         END DO
182      ELSE                                         ! Constant Volume case ==>> Concentration dilution effect
183         DO jj = 2, jpj
184            DO ji = fs_2, fs_jpim1   ! vector opt.
185               ! temperature : heat flux
186               sbc_tsc(ji,jj,jp_tem) = r1_rau0_rcp * qns(ji,jj)                          &   ! non solar heat flux
187                  &                  + r1_rau0     * emp(ji,jj)  * tsn(ji,jj,1,jp_tem)       ! concent./dilut. effect
188               ! salinity    : salt flux + concent./dilut. effect (both in sfx)
189               sbc_tsc(ji,jj,jp_sal) = r1_rau0  * (  sfx(ji,jj)                          &   ! salt flux (freezing/melting)
190                  &                                + emp(ji,jj) * tsn(ji,jj,1,jp_sal) )      ! concent./dilut. effect
191            END DO
192         END DO
193         IF( iom_use('emp_x_sst') )   CALL iom_put( "emp_x_sst", emp (:,:) * tsn(:,:,1,jp_tem) )   ! c/d term on sst
194         IF( iom_use('emp_x_sss') )   CALL iom_put( "emp_x_sss", emp (:,:) * tsn(:,:,1,jp_sal) )   ! c/d term on sss
195      ENDIF
196      ! Concentration dilution effect on (t,s) due to evapouration, precipitation and qns, but not river runoff 
197      DO jn = 1, jpts
198         DO jj = 2, jpj
199            DO ji = fs_2, fs_jpim1   ! vector opt.
200               z1_e3t = zfact / fse3t(ji,jj,1)
201               tsa(ji,jj,1,jn) = tsa(ji,jj,1,jn) + ( sbc_tsc_b(ji,jj,jn) + sbc_tsc(ji,jj,jn) ) * z1_e3t
202            END DO
203         END DO
204      END DO
205      !                                          Write in the ocean restart file
206      !                                          *******************************
207      IF( lrst_oce ) THEN
208         IF(lwp) WRITE(numout,*)
209         IF(lwp) WRITE(numout,*) 'sbc : ocean surface tracer content forcing fields written in ocean restart file ',   &
210            &                    'at it= ', kt,' date= ', ndastp
211         IF(lwp) WRITE(numout,*) '~~~~'
212         IF(nn_timing == 2)  CALL timing_start('iom_rstput')
213         CALL iom_rstput( kt, nitrst, numrow, 'sbc_hc_b', sbc_tsc(:,:,jp_tem) )
214         CALL iom_rstput( kt, nitrst, numrow, 'sbc_sc_b', sbc_tsc(:,:,jp_sal) )
215         IF(nn_timing == 2)  CALL timing_stop('iom_rstput')
216      ENDIF
217      !
218      !
219      !----------------------------------------
220      !       Ice Shelf effects (ISF)
221      !     tbl treated as in Losh (2008) JGR
222      !----------------------------------------
223      !
224      IF( nn_isf > 0 ) THEN
225         zfact = 0.5_wp
226         DO jj = 2, jpj
227            DO ji = fs_2, fs_jpim1
228         
229               ikt = misfkt(ji,jj)
230               ikb = misfkb(ji,jj)
231   
232               ! level fully include in the ice shelf boundary layer
233               ! sign - because fwf sign of evapo (rnf sign of precip)
234               DO jk = ikt, ikb - 1
235               ! compute trend
236                  tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem)                                                &
237                     &           + zfact * ( risf_tsc_b(ji,jj,jp_tem) + risf_tsc(ji,jj,jp_tem) )             &
238                     &           * r1_hisf_tbl(ji,jj)
239               END DO
240   
241               ! level partially include in ice shelf boundary layer
242               ! compute trend
243               tsa(ji,jj,ikb,jp_tem) = tsa(ji,jj,ikb,jp_tem)                                                 &
244                  &              + zfact * ( risf_tsc_b(ji,jj,jp_tem) + risf_tsc(ji,jj,jp_tem) )             &
245                  &              * r1_hisf_tbl(ji,jj) * ralpha(ji,jj)
246
247            END DO
248         END DO
249         IF( lrst_oce ) THEN
250            IF(lwp) WRITE(numout,*)
251            IF(lwp) WRITE(numout,*) 'sbc : isf surface tracer content forcing fields written in ocean restart file ',   &
252               &                    'at it= ', kt,' date= ', ndastp
253            IF(lwp) WRITE(numout,*) '~~~~'
254            IF(nn_timing == 2)  CALL timing_start('iom_rstput')
255            CALL iom_rstput( kt, nitrst, numrow, 'fwf_isf_b', fwfisf(:,:)          )
256            CALL iom_rstput( kt, nitrst, numrow, 'isf_hc_b' , risf_tsc(:,:,jp_tem) )
257            CALL iom_rstput( kt, nitrst, numrow, 'isf_sc_b' , risf_tsc(:,:,jp_sal) )
258            IF(nn_timing == 2)  CALL timing_stop('iom_rstput')
259         ENDIF
260      END IF
261      !
262      !----------------------------------------
263      !        River Runoff effects
264      !----------------------------------------
265      !
266      IF( ln_rnf ) THEN         ! input of heat and salt due to river runoff
267         zfact = 0.5_wp
268         DO jj = 2, jpj 
269            DO ji = fs_2, fs_jpim1
270               IF( rnf(ji,jj) /= 0._wp ) THEN
271                  zdep = zfact / h_rnf(ji,jj)
272                  DO jk = 1, nk_rnf(ji,jj)
273                                        tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem)   &
274                                          &               +  ( rnf_tsc_b(ji,jj,jp_tem) + rnf_tsc(ji,jj,jp_tem) ) * zdep
275                     IF( ln_rnf_sal )   tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal)   &
276                                          &               +  ( rnf_tsc_b(ji,jj,jp_sal) + rnf_tsc(ji,jj,jp_sal) ) * zdep 
277                  END DO
278               ENDIF
279            END DO 
280         END DO 
281      ENDIF
282
283      IF( iom_use('rnf_x_sst') )   CALL iom_put( "rnf_x_sst", rnf*tsn(:,:,1,jp_tem) )   ! runoff term on sst
284      IF( iom_use('rnf_x_sss') )   CALL iom_put( "rnf_x_sss", rnf*tsn(:,:,1,jp_sal) )   ! runoff term on sss
285
286#if defined key_asminc
287! WARNING: THIS MAY WELL NOT BE REQUIRED - WE DON'T WANT TO CHANGE T&S BUT THIS MAY COMPENSATE ANOTHER TERM...
288! Rate of change in e3t for each level is ssh_iau*e3t_0/ht_0
289! Contribution to tsa should be rate of change in level / per m of ocean? (hence the division by fse3t_n)
290      IF( ln_sshinc ) THEN         ! input of heat and salt due to assimilation
291         DO jj = 2, jpj 
292            DO ji = fs_2, fs_jpim1
293               zdep = ssh_iau(ji,jj) / ( ht_0(ji,jj) + 1.0 - ssmask(ji, jj) )
294               DO jk = 1, jpkm1
295                  tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem)   &
296                                        &            + tsn(ji,jj,jk,jp_tem) * zdep * ( e3t_0(ji,jj,jk) / fse3t_n(ji,jj,jk) )
297                  tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal)   &
298                                        &            + tsn(ji,jj,jk,jp_sal) * zdep * ( e3t_0(ji,jj,jk) / fse3t_n(ji,jj,jk) )
299               END DO
300            END DO 
301         END DO 
302      ENDIF
303#endif
304 
305      IF( l_trdtra )   THEN                      ! send trends for further diagnostics
306         ztrdt(:,:,:) = tsa(:,:,:,jp_tem) - ztrdt(:,:,:)
307         ztrds(:,:,:) = tsa(:,:,:,jp_sal) - ztrds(:,:,:)
308         CALL trd_tra( kt, 'TRA', jp_tem, jptra_nsr, ztrdt )
309         CALL trd_tra( kt, 'TRA', jp_sal, jptra_nsr, ztrds )
310         CALL wrk_dealloc( jpi, jpj, jpk, ztrdt, ztrds ) 
311      ENDIF
312      !
313      IF(ln_ctl)   CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' sbc  - Ta: ', mask1=tmask,   &
314         &                       tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
315      !
316      IF( nn_timing == 1 )  CALL timing_stop('tra_sbc')
317      !
318   END SUBROUTINE tra_sbc
319
320   !!======================================================================
321END MODULE trasbc
Note: See TracBrowser for help on using the repository browser.