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/2015/dev_r5151_UKMO_ISF/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/2015/dev_r5151_UKMO_ISF/NEMOGCM/NEMO/OPA_SRC/TRA/trasbc.F90 @ 5624

Last change on this file since 5624 was 5624, checked in by mathiot, 9 years ago

UKMO_ISF : fix conservation issue based on the work of Jerome on runoff, simplification of trasbc (isf part only) and remove option to apply isf melting as volume flux or not

  • Property svn:keywords set to Id
File size: 15.2 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
36   IMPLICIT NONE
37   PRIVATE
38
39   PUBLIC   tra_sbc    ! routine called by step.F90
40
41   !! * Substitutions
42#  include "domzgr_substitute.h90"
43#  include "vectopt_loop_substitute.h90"
44   !!----------------------------------------------------------------------
45   !! NEMO/OPA 3.7 , NEMO Consortium (2014)
46   !! $Id$
47   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
48   !!----------------------------------------------------------------------
49CONTAINS
50
51   SUBROUTINE tra_sbc ( kt )
52      !!----------------------------------------------------------------------
53      !!                  ***  ROUTINE tra_sbc  ***
54      !!                   
55      !! ** Purpose :   Compute the tracer surface boundary condition trend of
56      !!      (flux through the interface, concentration/dilution effect)
57      !!      and add it to the general trend of tracer equations.
58      !!
59      !! ** Method :
60      !!      Following Roullet and Madec (2000), the air-sea flux can be divided
61      !!      into three effects: (1) Fext, external forcing;
62      !!      (2) Fwi, concentration/dilution effect due to water exchanged
63      !!         at the surface by evaporation, precipitations and runoff (E-P-R);
64      !!      (3) Fwe, tracer carried with the water that is exchanged.
65      !!            - salinity    : salt flux only due to freezing/melting
66      !!            sa = sa +  sfx / rau0 / e3t  for k=1
67      !!
68      !!      Fext, flux through the air-sea interface for temperature and salt:
69      !!            - temperature : heat flux q (w/m2). If penetrative solar
70      !!         radiation q is only the non solar part of the heat flux, the
71      !!         solar part is added in traqsr.F routine.
72      !!            ta = ta + q /(rau0 rcp e3t)  for k=1
73      !!            - salinity    : no salt flux
74      !!
75      !!      The formulation for Fwb and Fwi vary according to the free
76      !!      surface formulation (linear or variable volume).
77      !!      * Linear free surface
78      !!            The surface freshwater flux modifies the ocean volume
79      !!         and thus the concentration of a tracer and the temperature.
80      !!         First order of the effect of surface freshwater exchange
81      !!         for salinity, it can be neglected on temperature (especially
82      !!         as the temperature of precipitations and runoffs is usually
83      !!         unknown).
84      !!            - temperature : we assume that the temperature of both
85      !!         precipitations and runoffs is equal to the SST, thus there
86      !!         is no additional flux since in this case, the concentration
87      !!         dilution effect is balanced by the net heat flux associated
88      !!         to the freshwater exchange (Fwe+Fwi=0):
89      !!            (Tp P - Te E) + SST (P-E) = 0 when Tp=Te=SST
90      !!            - salinity    : evaporation, precipitation and runoff
91      !!         water has a zero salinity  but there is a salt flux due to
92      !!         freezing/melting, thus:
93      !!            sa = sa + emp * sn / rau0 / e3t   for k=1
94      !!                    + sfx    / rau0 / e3t
95      !!         where emp, the surface freshwater budget (evaporation minus
96      !!         precipitation minus runoff) given in kg/m2/s is divided
97      !!         by rau0 (density of sea water) to obtain m/s.   
98      !!         Note: even though Fwe does not appear explicitly for
99      !!         temperature in this routine, the heat carried by the water
100      !!         exchanged through the surface is part of the total heat flux
101      !!         forcing and must be taken into account in the global heat
102      !!         balance).
103      !!      * nonlinear free surface (variable volume, lk_vvl)
104      !!         contrary to the linear free surface case, Fwi is properly
105      !!         taken into account by using the true layer thicknesses to       
106      !!         calculate tracer content and advection. There is no need to
107      !!         deal with it in this routine.
108      !!           - temperature: Fwe=SST (P-E+R) is added to Fext.
109      !!           - salinity:  Fwe = 0, there is no surface flux of salt.
110      !!
111      !! ** Action  : - Update the 1st level of (ta,sa) with the trend associated
112      !!                with the tracer surface boundary condition
113      !!              - send trends to trdtra module (l_trdtra=T)
114      !!----------------------------------------------------------------------
115      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
116      !!
117      INTEGER  ::   ji, jj, jk, jn           ! dummy loop indices 
118      INTEGER  ::   ikt, ikb 
119      REAL(wp) ::   zfact, z1_e3t, zdep
120      REAL(wp) ::   zt_frz, zpress
121      REAL(wp), POINTER, DIMENSION(:,:,:) ::  ztrdt, ztrds
122      !!----------------------------------------------------------------------
123      !
124      IF( nn_timing == 1 )  CALL timing_start('tra_sbc')
125      !
126      IF( kt == nit000 ) THEN
127         IF(lwp) WRITE(numout,*)
128         IF(lwp) WRITE(numout,*) 'tra_sbc : TRAcer Surface Boundary Condition'
129         IF(lwp) WRITE(numout,*) '~~~~~~~ '
130      ENDIF
131
132      IF( l_trdtra ) THEN                    !* Save ta and sa trends
133         CALL wrk_alloc( jpi, jpj, jpk, ztrdt, ztrds ) 
134         ztrdt(:,:,:) = tsa(:,:,:,jp_tem)
135         ztrds(:,:,:) = tsa(:,:,:,jp_sal)
136      ENDIF
137
138!!gm      IF( .NOT.ln_traqsr )   qsr(:,:) = 0.e0   ! no solar radiation penetration
139      IF( .NOT.ln_traqsr ) THEN     ! no solar radiation penetration
140         qns(:,:) = qns(:,:) + qsr(:,:)      ! total heat flux in qns
141         qsr(:,:) = 0.e0                     ! qsr set to zero
142      ENDIF
143
144      !----------------------------------------
145      !        EMP, SFX and QNS effects
146      !----------------------------------------
147      !                                          Set before sbc tracer content fields
148      !                                          ************************************
149      IF( kt == nit000 ) THEN                      ! Set the forcing field at nit000 - 1
150         !                                         ! -----------------------------------
151         IF( ln_rstart .AND.    &                     ! Restart: read in restart file
152              & iom_varid( numror, 'sbc_hc_b', ldstop = .FALSE. ) > 0 ) THEN
153            IF(lwp) WRITE(numout,*) '          nit000-1 surface tracer content forcing fields red in the restart file'
154            zfact = 0.5_wp
155            CALL iom_get( numror, jpdom_autoglo, 'sbc_hc_b', sbc_tsc_b(:,:,jp_tem) )   ! before heat content sbc trend
156            CALL iom_get( numror, jpdom_autoglo, 'sbc_sc_b', sbc_tsc_b(:,:,jp_sal) )   ! before salt content sbc trend
157         ELSE                                         ! No restart or restart not found: Euler forward time stepping
158            zfact = 1._wp
159            sbc_tsc_b(:,:,:) = 0._wp
160         ENDIF
161      ELSE                                         ! Swap of forcing fields
162         !                                         ! ----------------------
163         zfact = 0.5_wp
164         sbc_tsc_b(:,:,:) = sbc_tsc(:,:,:)
165      ENDIF
166      !                                          Compute now sbc tracer content fields
167      !                                          *************************************
168
169                                                   ! Concentration dilution effect on (t,s) due to 
170                                                   ! evaporation, precipitation and qns, but not river runoff
171                                               
172      IF( lk_vvl ) THEN                            ! Variable Volume case  ==>> heat content of mass flux is in qns
173         DO jj = 1, jpj
174            DO ji = 1, jpi 
175               sbc_tsc(ji,jj,jp_tem) = r1_rau0_rcp * qns(ji,jj)                              ! non solar heat flux
176               sbc_tsc(ji,jj,jp_sal) = r1_rau0     * sfx(ji,jj)                              ! salt flux due to freezing/melting
177            END DO
178         END DO
179      ELSE                                         ! Constant Volume case ==>> Concentration dilution effect
180         DO jj = 2, jpj
181            DO ji = fs_2, fs_jpim1   ! vector opt.
182               ! temperature : heat flux
183               sbc_tsc(ji,jj,jp_tem) = r1_rau0_rcp * qns(ji,jj)                          &   ! non solar heat flux
184                  &                  + r1_rau0     * emp(ji,jj)  * tsn(ji,jj,1,jp_tem)       ! concent./dilut. effect
185               ! salinity    : salt flux + concent./dilut. effect (both in sfx)
186               sbc_tsc(ji,jj,jp_sal) = r1_rau0  * (  sfx(ji,jj)                          &   ! salt flux (freezing/melting)
187                  &                                + emp(ji,jj) * tsn(ji,jj,1,jp_sal) )      ! concent./dilut. effect
188            END DO
189         END DO
190         IF( iom_use('emp_x_sst') )   CALL iom_put( "emp_x_sst", emp (:,:) * tsn(:,:,1,jp_tem) )   ! c/d term on sst
191         IF( iom_use('emp_x_sss') )   CALL iom_put( "emp_x_sss", emp (:,:) * tsn(:,:,1,jp_sal) )   ! c/d term on sss
192      ENDIF
193      ! Concentration dilution effect on (t,s) due to evapouration, precipitation and qns, but not river runoff 
194      DO jn = 1, jpts
195         DO jj = 2, jpj
196            DO ji = fs_2, fs_jpim1   ! vector opt.
197               z1_e3t = zfact / fse3t(ji,jj,1)
198               tsa(ji,jj,1,jn) = tsa(ji,jj,1,jn) + ( sbc_tsc_b(ji,jj,jn) + sbc_tsc(ji,jj,jn) ) * z1_e3t
199            END DO
200         END DO
201      END DO
202      !                                          Write in the ocean restart file
203      !                                          *******************************
204      IF( lrst_oce ) THEN
205         IF(lwp) WRITE(numout,*)
206         IF(lwp) WRITE(numout,*) 'sbc : ocean surface tracer content forcing fields written in ocean restart file ',   &
207            &                    'at it= ', kt,' date= ', ndastp
208         IF(lwp) WRITE(numout,*) '~~~~'
209         CALL iom_rstput( kt, nitrst, numrow, 'sbc_hc_b', sbc_tsc(:,:,jp_tem) )
210         CALL iom_rstput( kt, nitrst, numrow, 'sbc_sc_b', sbc_tsc(:,:,jp_sal) )
211      ENDIF
212      !
213      !
214      !----------------------------------------
215      !       Ice Shelf effects (ISF)
216      !     tbl treated as in Losh (2008) JGR
217      !----------------------------------------
218      !
219      IF( ln_isf ) THEN
220         zfact = 0.5_wp
221         DO jj = 2, jpj
222            DO ji = fs_2, fs_jpim1
223         
224               ikt = misfkt(ji,jj)
225               ikb = misfkb(ji,jj)
226   
227               ! level fully include in the ice shelf boundary layer
228               ! sign - because fwf sign of evapo (rnf sign of precip)
229               DO jk = ikt, ikb - 1
230               ! compute trend
231                  tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem)                                                &
232                     &           + zfact * ( risf_tsc_b(ji,jj,jp_tem) + risf_tsc(ji,jj,jp_tem) )             &
233                     &           * r1_hisf_tbl(ji,jj)
234               END DO
235   
236               ! level partially include in ice shelf boundary layer
237               ! compute trend
238               tsa(ji,jj,ikb,jp_tem) = tsa(ji,jj,ikb,jp_tem)                                                 &
239                  &              + zfact * ( risf_tsc_b(ji,jj,jp_tem) + risf_tsc(ji,jj,jp_tem) )             &
240                  &              * r1_hisf_tbl(ji,jj) * ralpha(ji,jj)
241
242            END DO
243         END DO
244         IF( lrst_oce ) THEN
245            IF(lwp) WRITE(numout,*)
246            IF(lwp) WRITE(numout,*) 'sbc : isf surface tracer content forcing fields written in ocean restart file ',   &
247               &                    'at it= ', kt,' date= ', ndastp
248            IF(lwp) WRITE(numout,*) '~~~~'
249            CALL iom_rstput( kt, nitrst, numrow, 'fwf_isf_b', fwfisf(:,:)          )
250            CALL iom_rstput( kt, nitrst, numrow, 'isf_hc_b' , risf_tsc(:,:,jp_tem) )
251            CALL iom_rstput( kt, nitrst, numrow, 'isf_sc_b' , risf_tsc(:,:,jp_sal) )
252         ENDIF
253      END IF
254      !
255      !----------------------------------------
256      !        River Runoff effects
257      !----------------------------------------
258      !
259      IF( ln_rnf ) THEN         ! input of heat and salt due to river runoff
260         zfact = 0.5_wp
261         DO jj = 2, jpj 
262            DO ji = fs_2, fs_jpim1
263               IF( rnf(ji,jj) /= 0._wp ) THEN
264                  zdep = zfact / h_rnf(ji,jj)
265                  DO jk = 1, nk_rnf(ji,jj)
266                                        tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem)   &
267                                          &               +  ( rnf_tsc_b(ji,jj,jp_tem) + rnf_tsc(ji,jj,jp_tem) ) * zdep
268                     IF( ln_rnf_sal )   tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal)   &
269                                          &               +  ( rnf_tsc_b(ji,jj,jp_sal) + rnf_tsc(ji,jj,jp_sal) ) * zdep 
270                  END DO
271               ENDIF
272            END DO 
273         END DO 
274      ENDIF
275 
276      IF( l_trdtra )   THEN                      ! send trends for further diagnostics
277         ztrdt(:,:,:) = tsa(:,:,:,jp_tem) - ztrdt(:,:,:)
278         ztrds(:,:,:) = tsa(:,:,:,jp_sal) - ztrds(:,:,:)
279         CALL trd_tra( kt, 'TRA', jp_tem, jptra_nsr, ztrdt )
280         CALL trd_tra( kt, 'TRA', jp_sal, jptra_nsr, ztrds )
281         CALL wrk_dealloc( jpi, jpj, jpk, ztrdt, ztrds ) 
282      ENDIF
283      !
284      IF(ln_ctl)   CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' sbc  - Ta: ', mask1=tmask,   &
285         &                       tab3d_2=tsa(:,:,:,jp_sal), clinfo2=       ' Sa: ', mask2=tmask, clinfo3='tra' )
286      !
287      IF( nn_timing == 1 )  CALL timing_stop('tra_sbc')
288      !
289   END SUBROUTINE tra_sbc
290
291   !!======================================================================
292END MODULE trasbc
Note: See TracBrowser for help on using the repository browser.