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.
trcbc.F90 in NEMO/branches/2019/dev_r12072_MERGE_OPTION2_2019/src/TOP – NEMO

source: NEMO/branches/2019/dev_r12072_MERGE_OPTION2_2019/src/TOP/trcbc.F90 @ 12708

Last change on this file since 12708 was 12327, checked in by cetlod, 4 years ago

dev_r12072_MERGE_OPTION2_2019: minor fix to avoid compilation error when using an OCE+TOP configuration without key_top enabled

  • Property svn:keywords set to Id
File size: 23.2 KB
RevLine 
[4058]1MODULE trcbc
2   !!======================================================================
[6140]3   !!                     ***  MODULE  trcbc  ***
[4058]4   !! TOP :  module for passive tracer boundary conditions
5   !!=====================================================================
[7646]6   !! History :  3.5 !  2014 (M. Vichi, T. Lovato)  Original
7   !!            3.6 !  2015 (T . Lovato) Revision and BDY support
8   !!            4.0 !  2016 (T . Lovato) Include application of sbc and cbc
[4058]9   !!----------------------------------------------------------------------
[7646]10   !!   trc_bc       :  Apply tracer Boundary Conditions
[4058]11   !!----------------------------------------------------------------------
12   USE par_trc       !  passive tracers parameters
13   USE oce_trc       !  shared variables between ocean and passive tracers
14   USE trc           !  passive tracers common variables
15   USE iom           !  I/O manager
16   USE lib_mpp       !  MPP library
17   USE fldread       !  read input fields
[7646]18   USE bdy_oce,  ONLY: ln_bdy, nb_bdy , idx_bdy, ln_coords_file, rn_time_dmp, rn_time_dmp_out
[4058]19
20   IMPLICIT NONE
21   PRIVATE
22
[7646]23   PUBLIC   trc_bc         ! called in trcstp.F90 or within TOP modules
24   PUBLIC   trc_bc_ini     ! called in trcini.F90
[4058]25
[6140]26   INTEGER  , SAVE, PUBLIC                             :: nb_trcobc    ! number of tracers with open BC
27   INTEGER  , SAVE, PUBLIC                             :: nb_trcsbc    ! number of tracers with surface BC
28   INTEGER  , SAVE, PUBLIC                             :: nb_trccbc    ! number of tracers with coastal BC
[4058]29   INTEGER  , SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: n_trc_indobc ! index of tracer with OBC data
30   INTEGER  , SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: n_trc_indsbc ! index of tracer with SBC data
31   INTEGER  , SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: n_trc_indcbc ! index of tracer with CBC data
[6140]32   REAL(wp) , SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: rf_trsfac    ! multiplicative factor for SBC tracer values
33   TYPE(FLD), SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: sf_trcsbc    ! structure of data input SBC (file informations, fields read)
34   REAL(wp) , SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: rf_trcfac    ! multiplicative factor for CBC tracer values
35   TYPE(FLD), SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: sf_trccbc    ! structure of data input CBC (file informations, fields read)
36   REAL(wp) , SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: rf_trofac    ! multiplicative factor for OBCtracer values
[9800]37#if defined key_agrif
[9788]38   TYPE(FLD), SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: sf_trcobc    ! structure of data input OBC (file informations, fields read)
[9800]39#else
40   TYPE(FLD), SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:), TARGET  :: sf_trcobc
41#endif
[4058]42
[12327]43#if defined key_top
44   !!----------------------------------------------------------------------
45   !!   'key_top'                                                TOP model
46   !!----------------------------------------------------------------------
47
[7646]48   !! * Substitutions
49#  include "vectopt_loop_substitute.h90"
[4058]50   !!----------------------------------------------------------------------
[9598]51   !! NEMO/TOP 4.0 , NEMO Consortium (2018)
[6140]52   !! $Id$
[10068]53   !! Software governed by the CeCILL license (see ./LICENSE)
[4058]54   !!----------------------------------------------------------------------
55CONTAINS
56
[7646]57   SUBROUTINE trc_bc_ini( ntrc )
[4058]58      !!----------------------------------------------------------------------
[7646]59      !!                   ***  ROUTINE trc_bc_ini  ***
[4058]60      !!                   
61      !! ** Purpose :   initialisation of passive tracer BC data
62      !!
63      !! ** Method  : - Read namtsd namelist
64      !!              - allocates passive tracer BC data structure
65      !!----------------------------------------------------------------------
[9169]66      INTEGER,INTENT(in) :: ntrc                           ! number of tracers
67      !
[6140]68      INTEGER            :: jl, jn , ib, ibd, ii, ij, ik   ! dummy loop indices
[4058]69      INTEGER            :: ierr0, ierr1, ierr2, ierr3     ! temporary integers
[6140]70      INTEGER            :: ios                            ! Local integer output status for namelist read
71      INTEGER            :: nblen, igrd                    ! support arrays for BDY
[4058]72      CHARACTER(len=100) :: clndta, clntrc
73      !
[6140]74      CHARACTER(len=100) :: cn_dir_sbc, cn_dir_cbc, cn_dir_obc
[4058]75      TYPE(FLD_N), ALLOCATABLE, DIMENSION(:) :: slf_i  ! local array of namelist informations on the fields to read
76      TYPE(FLD_N), DIMENSION(jpmaxtrc) :: sn_trcobc    ! open
77      TYPE(FLD_N), DIMENSION(jpmaxtrc) :: sn_trcsbc    ! surface
78      TYPE(FLD_N), DIMENSION(jpmaxtrc) :: sn_trccbc    ! coastal
79      REAL(wp)   , DIMENSION(jpmaxtrc) :: rn_trofac    ! multiplicative factor for tracer values
80      REAL(wp)   , DIMENSION(jpmaxtrc) :: rn_trsfac    ! multiplicative factor for tracer values
81      REAL(wp)   , DIMENSION(jpmaxtrc) :: rn_trcfac    ! multiplicative factor for tracer values
82      !!
[7646]83      NAMELIST/namtrc_bc/ cn_dir_obc, sn_trcobc, rn_trofac, cn_dir_sbc, sn_trcsbc, rn_trsfac, & 
[12110]84                        & cn_dir_cbc, sn_trccbc, rn_trcfac, ln_rnf_ctl, rn_sbc_time, rn_cbc_time
[6140]85      NAMELIST/namtrc_bdy/ cn_trc_dflt, cn_trc, nn_trcdmp_bdy
[4058]86      !!----------------------------------------------------------------------
87      !
[6140]88      IF( lwp ) THEN
[9169]89         WRITE(numout,*)
[7646]90         WRITE(numout,*) 'trc_bc_ini : Tracers Boundary Conditions (BC)'
[6140]91         WRITE(numout,*) '~~~~~~~~~~~ '
92      ENDIF
[4058]93      !  Initialisation and local array allocation
[9169]94      ierr0 = 0   ;   ierr1 = 0   ;   ierr2 = 0   ;   ierr3 = 0 
[4058]95      ALLOCATE( slf_i(ntrc), STAT=ierr0 )
96      IF( ierr0 > 0 ) THEN
[7646]97         CALL ctl_stop( 'trc_bc_ini: unable to allocate local slf_i' )   ;   RETURN
[4058]98      ENDIF
99
100      ! Compute the number of tracers to be initialised with open, surface and boundary data
101      ALLOCATE( n_trc_indobc(ntrc), STAT=ierr0 )
102      IF( ierr0 > 0 ) THEN
[7646]103         CALL ctl_stop( 'trc_bc_ini: unable to allocate n_trc_indobc' )   ;   RETURN
[4058]104      ENDIF
[9169]105      nb_trcobc       = 0
[4058]106      n_trc_indobc(:) = 0
107      !
108      ALLOCATE( n_trc_indsbc(ntrc), STAT=ierr0 )
109      IF( ierr0 > 0 ) THEN
[7646]110         CALL ctl_stop( 'trc_bc_ini: unable to allocate n_trc_indsbc' )   ;   RETURN
[4058]111      ENDIF
[9169]112      nb_trcsbc       = 0
[4058]113      n_trc_indsbc(:) = 0
114      !
115      ALLOCATE( n_trc_indcbc(ntrc), STAT=ierr0 )
116      IF( ierr0 > 0 ) THEN
[7646]117         CALL ctl_stop( 'trc_bc_ini: unable to allocate n_trc_indcbc' )   ;   RETURN
[4058]118      ENDIF
[9169]119      nb_trccbc       = 0
[4058]120      n_trc_indcbc(:) = 0
121      !
[6140]122      ! Read Boundary Conditions Namelists
123      READ  ( numnat_ref, namtrc_bc, IOSTAT = ios, ERR = 901)
[11536]124901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namtrc_bc in reference namelist' )
[6140]125      READ  ( numnat_cfg, namtrc_bc, IOSTAT = ios, ERR = 902 )
[11536]126902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namtrc_bc in configuration namelist' )
[6140]127      IF(lwm) WRITE ( numont, namtrc_bc )
128
[7646]129      IF ( ln_bdy ) THEN
130         READ  ( numnat_ref, namtrc_bdy, IOSTAT = ios, ERR = 903)
[11536]131903      IF( ios /= 0 )   CALL ctl_nam ( ios , 'namtrc_bdy in reference namelist' )
132         ! make sur that all elements of the namelist variables have a default definition from namelist_ref
133         cn_trc     (2:jp_bdy) = cn_trc     (1)
134         cn_trc_dflt(2:jp_bdy) = cn_trc_dflt(1)
[7646]135         READ  ( numnat_cfg, namtrc_bdy, IOSTAT = ios, ERR = 904 )
[11536]136904      IF( ios >  0 )   CALL ctl_nam ( ios , 'namtrc_bdy in configuration namelist' )
[7646]137         IF(lwm) WRITE ( numont, namtrc_bdy )
138     
139         ! setup up preliminary informations for BDY structure
140         DO jn = 1, ntrc
141            DO ib = 1, nb_bdy
142               ! Set type of obc in BDY data structure (around here we may plug user override of obc type from nml)
[9169]143               IF ( ln_trc_obc(jn) ) THEN   ;   trcdta_bdy(jn,ib)%cn_obc = TRIM( cn_trc     (ib) )
144               ELSE                         ;   trcdta_bdy(jn,ib)%cn_obc = TRIM( cn_trc_dflt(ib) )
[7646]145               ENDIF
146               ! set damping use in BDY data structure
147               trcdta_bdy(jn,ib)%dmp = .false.
[9169]148               IF(nn_trcdmp_bdy(ib) == 1 .AND. ln_trc_obc(jn) )   trcdta_bdy(jn,ib)%dmp = .true.
149               IF(nn_trcdmp_bdy(ib) == 2                      )   trcdta_bdy(jn,ib)%dmp = .true.
150               IF(trcdta_bdy(jn,ib)%cn_obc == 'frs' .AND. nn_trcdmp_bdy(ib) /= 0 )  &
151                   & CALL ctl_stop( 'trc_bc_ini: Use FRS OR relaxation' )
152               IF(  .NOT.( 0 < nn_trcdmp_bdy(ib)  .AND.  nn_trcdmp_bdy(ib) <= 2 )  )   &
153                   & CALL ctl_stop( 'trc_bc_ini: Not a valid option for nn_trcdmp_bdy. Allowed: 0,1,2.' )
154            END DO
155         END DO
[7646]156      ELSE
157         ! Force all tracers OBC to false if bdy not used
158         ln_trc_obc = .false.
159      ENDIF
[6140]160
161      ! compose BC data indexes
162      DO jn = 1, ntrc
[4058]163         IF( ln_trc_obc(jn) ) THEN
[9169]164             nb_trcobc       = nb_trcobc + 1   ;   n_trc_indobc(jn) = nb_trcobc
[4058]165         ENDIF
166         IF( ln_trc_sbc(jn) ) THEN
[9169]167             nb_trcsbc       = nb_trcsbc + 1   ;   n_trc_indsbc(jn) = nb_trcsbc
[4058]168         ENDIF
169         IF( ln_trc_cbc(jn) ) THEN
[9169]170             nb_trccbc       = nb_trccbc + 1   ;   n_trc_indcbc(jn) = nb_trccbc
[4058]171         ENDIF
[9169]172      END DO
[4058]173
[6140]174      ! Print summmary of Boundary Conditions
175      IF( lwp ) THEN
[9169]176         WRITE(numout,*)
[6140]177         WRITE(numout,'(a,i3)') '   Total tracers to be initialized with SURFACE BCs data:', nb_trcsbc
178         IF ( nb_trcsbc > 0 ) THEN
179            WRITE(numout,*) '   #trc        NAME        Boundary     Mult.Fact. '
180            DO jn = 1, ntrc
181               IF ( ln_trc_sbc(jn) ) WRITE(numout,9001) jn, TRIM( sn_trcsbc(jn)%clvar ), 'SBC', rn_trsfac(jn)
[9169]182            END DO
[6140]183         ENDIF
184         WRITE(numout,'(2a)') '   SURFACE BC data repository : ', TRIM(cn_dir_sbc)
[9169]185         !
186         WRITE(numout,*)
[6140]187         WRITE(numout,'(a,i3)') '   Total tracers to be initialized with COASTAL BCs data:', nb_trccbc
[7646]188         IF( nb_trccbc > 0 ) THEN
[6140]189            WRITE(numout,*) '   #trc        NAME        Boundary     Mult.Fact. '
190            DO jn = 1, ntrc
191               IF ( ln_trc_cbc(jn) ) WRITE(numout, 9001) jn, TRIM( sn_trccbc(jn)%clvar ), 'CBC', rn_trcfac(jn)
[9169]192            END DO
[6140]193         ENDIF
194         WRITE(numout,'(2a)') '   COASTAL BC data repository : ', TRIM(cn_dir_cbc)
[9169]195         IF( .NOT.ln_rnf .OR. .NOT.ln_linssh )   ln_rnf_ctl = .FALSE.
[9119]196         IF( ln_rnf_ctl )  WRITE(numout,'(a)') &
197              &            ' -> Remove runoff dilution effect on tracers with absent river load (ln_rnf_ctl = .TRUE.)'
[9169]198         WRITE(numout,*)
[6140]199         WRITE(numout,'(a,i3)') '   Total tracers to be initialized with OPEN BCs data:', nb_trcobc
[7646]200
201         IF( ln_bdy .AND. nb_trcobc > 0 ) THEN
[6140]202            WRITE(numout,*) '   #trc        NAME        Boundary     Mult.Fact.   OBC Settings'
203            DO jn = 1, ntrc
[9119]204               IF (       ln_trc_obc(jn) )  WRITE(numout, 9001) jn, TRIM( sn_trcobc(jn)%clvar ), 'OBC', rn_trofac(jn), &
205                    &                                           (trcdta_bdy(jn,ib)%cn_obc,ib=1,nb_bdy)
206               IF ( .NOT. ln_trc_obc(jn) )  WRITE(numout, 9002) jn, 'Set data to IC and use default condition'       , &
207                    &                                           (trcdta_bdy(jn,ib)%cn_obc,ib=1,nb_bdy)
[9169]208            END DO
[6140]209            WRITE(numout,*) ' '
210            DO ib = 1, nb_bdy
[9169]211               IF(nn_trcdmp_bdy(ib) == 0) WRITE(numout,9003) '   Boundary ', ib, &
212                  &                                          ' -> NO damping of tracers'
213               IF(nn_trcdmp_bdy(ib) == 1) WRITE(numout,9003) '   Boundary ', ib, &
214                  &                                          ' -> damping ONLY for tracers with external data provided'
215               IF(nn_trcdmp_bdy(ib) == 2) WRITE(numout,9003) '   Boundary ', ib, &
216                  &                                          ' -> damping of ALL tracers'
217               IF(nn_trcdmp_bdy(ib) >  0) THEN
[6140]218                   WRITE(numout,9003) '     USE damping parameters from nambdy for boundary ', ib,' : '
[9169]219                   WRITE(numout,'(a,f10.2,a)') '     - Inflow damping time scale  : ',rn_time_dmp    (ib),' days'
[6140]220                   WRITE(numout,'(a,f10.2,a)') '     - Outflow damping time scale : ',rn_time_dmp_out(ib),' days'
[9169]221               ENDIF
222            END DO
[6140]223         ENDIF
[9169]224         !
[6140]225         WRITE(numout,'(2a)') '   OPEN BC data repository : ', TRIM(cn_dir_obc)
[4058]226      ENDIF
[6140]2279001  FORMAT(2x,i5, 3x, a15, 3x, a5, 6x, e11.3, 4x, 10a13)
2289002  FORMAT(2x,i5, 3x, a41, 3x, 10a13)
2299003  FORMAT(a, i5, a)
[4058]230      !
[9169]231      !
[4058]232      ! OPEN Lateral boundary conditions
[7646]233      IF( ln_bdy .AND. nb_trcobc > 0 ) THEN
[11536]234         ALLOCATE ( sf_trcobc(nb_trcobc), rf_trofac(nb_trcobc), STAT=ierr1 )
[4058]235         IF( ierr1 > 0 ) THEN
[7646]236            CALL ctl_stop( 'trc_bc_ini: unable to allocate sf_trcobc structure' )   ;   RETURN
[4058]237         ENDIF
[9169]238         !
[6140]239         igrd = 1                       ! Everything is at T-points here
[9169]240         !
[4058]241         DO jn = 1, ntrc
[6140]242            DO ib = 1, nb_bdy
[9169]243               !
[6140]244               nblen = idx_bdy(ib)%nblen(igrd)
[9169]245               !
246               IF( ln_trc_obc(jn) ) THEN     !* Initialise from external data *!
[6140]247                  jl = n_trc_indobc(jn)
248                  slf_i(jl)    = sn_trcobc(jn)
249                  rf_trofac(jl) = rn_trofac(jn)
[9169]250                                                ALLOCATE( sf_trcobc(jl)%fnow(nblen,1,jpk)   , STAT=ierr2 )
251                  IF( sn_trcobc(jn)%ln_tint )   ALLOCATE( sf_trcobc(jl)%fdta(nblen,1,jpk,2) , STAT=ierr3 )
[6140]252                  IF( ierr2 + ierr3 > 0 ) THEN
[7646]253                    CALL ctl_stop( 'trc_bc_ini : unable to allocate passive tracer OBC data arrays' )   ;   RETURN
[6140]254                  ENDIF
255                  trcdta_bdy(jn,ib)%trc => sf_trcobc(jl)%fnow(:,1,:)
256                  trcdta_bdy(jn,ib)%rn_fac = rf_trofac(jl)
[9169]257               ELSE                          !* Initialise obc arrays from initial conditions *!
[6140]258                  ALLOCATE ( trcdta_bdy(jn,ib)%trc(nblen,jpk) )
259                  DO ibd = 1, nblen
260                     DO ik = 1, jpkm1
261                        ii = idx_bdy(ib)%nbi(ibd,igrd)
262                        ij = idx_bdy(ib)%nbj(ibd,igrd)
263                        trcdta_bdy(jn,ib)%trc(ibd,ik) = trn(ii,ij,ik,jn) * tmask(ii,ij,ik)
264                     END DO
265                  END DO
266                  trcdta_bdy(jn,ib)%rn_fac = 1._wp
[4058]267               ENDIF
[9169]268            END DO
269         END DO
270         !
[7646]271         CALL fld_fill( sf_trcobc, slf_i, cn_dir_obc, 'trc_bc_ini', 'Passive tracer OBC data', 'namtrc_bc' )
[11536]272         DO jn = 1, ntrc   ! define imap pointer, must be done after the call to fld_fill
273            DO ib = 1, nb_bdy
274               IF( ln_trc_obc(jn) ) THEN     !* Initialise from external data *!
275                  jl = n_trc_indobc(jn)
276                  sf_trcobc(jl)%imap => idx_bdy(ib)%nbmap(1:idx_bdy(ib)%nblen(igrd),igrd)
277               ENDIF
278            END DO
279         END DO
280         !
[4058]281      ENDIF
[7646]282
[4058]283      ! SURFACE Boundary conditions
284      IF( nb_trcsbc > 0 ) THEN       !  allocate only if the number of tracer to initialise is greater than zero
285         ALLOCATE( sf_trcsbc(nb_trcsbc), rf_trsfac(nb_trcsbc), STAT=ierr1 )
286         IF( ierr1 > 0 ) THEN
[7646]287            CALL ctl_stop( 'trc_bc_ini: unable to allocate  sf_trcsbc structure' )   ;   RETURN
[4058]288         ENDIF
289         !
290         DO jn = 1, ntrc
291            IF( ln_trc_sbc(jn) ) THEN      ! update passive tracers arrays with input data read from file
292               jl = n_trc_indsbc(jn)
293               slf_i(jl)    = sn_trcsbc(jn)
294               rf_trsfac(jl) = rn_trsfac(jn)
295                                            ALLOCATE( sf_trcsbc(jl)%fnow(jpi,jpj,1)   , STAT=ierr2 )
296               IF( sn_trcsbc(jn)%ln_tint )  ALLOCATE( sf_trcsbc(jl)%fdta(jpi,jpj,1,2) , STAT=ierr3 )
297               IF( ierr2 + ierr3 > 0 ) THEN
[7646]298                 CALL ctl_stop( 'trc_bc_ini : unable to allocate passive tracer SBC data arrays' )   ;   RETURN
[4058]299               ENDIF
300            ENDIF
301            !   
[9169]302         END DO
[4058]303         !                         ! fill sf_trcsbc with slf_i and control print
[7646]304         CALL fld_fill( sf_trcsbc, slf_i, cn_dir_sbc, 'trc_bc_ini', 'Passive tracer SBC data', 'namtrc_bc' )
[4058]305         !
306      ENDIF
307      !
308      ! COSTAL Boundary conditions
309      IF( nb_trccbc > 0 ) THEN       !  allocate only if the number of tracer to initialise is greater than zero
310         ALLOCATE( sf_trccbc(nb_trccbc), rf_trcfac(nb_trccbc), STAT=ierr1 )
311         IF( ierr1 > 0 ) THEN
312            CALL ctl_stop( 'trc_bc_ini: unable to allocate  sf_trccbc structure' )   ;   RETURN
313         ENDIF
314         !
315         DO jn = 1, ntrc
316            IF( ln_trc_cbc(jn) ) THEN      ! update passive tracers arrays with input data read from file
317               jl = n_trc_indcbc(jn)
318               slf_i(jl)    = sn_trccbc(jn)
319               rf_trcfac(jl) = rn_trcfac(jn)
320                                            ALLOCATE( sf_trccbc(jl)%fnow(jpi,jpj,1)   , STAT=ierr2 )
321               IF( sn_trccbc(jn)%ln_tint )  ALLOCATE( sf_trccbc(jl)%fdta(jpi,jpj,1,2) , STAT=ierr3 )
322               IF( ierr2 + ierr3 > 0 ) THEN
323                 CALL ctl_stop( 'trc_bc_ini : unable to allocate passive tracer CBC data arrays' )   ;   RETURN
324               ENDIF
325            ENDIF
326            !   
[9169]327         END DO
[4058]328         !                         ! fill sf_trccbc with slf_i and control print
[7646]329         CALL fld_fill( sf_trccbc, slf_i, cn_dir_cbc, 'trc_bc_ini', 'Passive tracer CBC data', 'namtrc_bc' )
[4058]330         !
331      ENDIF
[6140]332      !
[4058]333      DEALLOCATE( slf_i )          ! deallocate local field structure
[6140]334      !
[7646]335   END SUBROUTINE trc_bc_ini
[4058]336
337
[7646]338   SUBROUTINE trc_bc(kt, jit)
[4058]339      !!----------------------------------------------------------------------
[7646]340      !!                   ***  ROUTINE trc_bc  ***
[4058]341      !!
[7646]342      !! ** Purpose :  Apply Boundary Conditions data to tracers
[4058]343      !!
[7646]344      !! ** Method  :  1) Read BC inputs and update data structures using fldread
345      !!               2) Apply Boundary Conditions to tracers
[4058]346      !!----------------------------------------------------------------------
347      USE fldread
[9124]348      !!     
[9169]349      INTEGER, INTENT(in)           ::   kt    ! ocean time-step index
350      INTEGER, INTENT(in), OPTIONAL ::   jit   ! subcycle time-step index (for timesplitting option)
[7646]351      !!
352      INTEGER  :: ji, jj, jk, jn, jl             ! Loop index
353      REAL(wp) :: zfact, zrnf
[4058]354      !!---------------------------------------------------------------------
355      !
[9124]356      IF( ln_timing )   CALL timing_start('trc_bc')
[4058]357
[6140]358      IF( kt == nit000 .AND. lwp) THEN
359         WRITE(numout,*)
[7646]360         WRITE(numout,*) 'trc_bc : Surface boundary conditions for passive tracers.'
[9169]361         WRITE(numout,*) '~~~~~~~ '
[4058]362      ENDIF
363
[7646]364      ! 1. Update Boundary conditions data
[9169]365      IF( PRESENT(jit) ) THEN 
366         !
[12256]367         ! BDY: use pt_offset=0.5 as applied at the end of the step and fldread is referenced at the middle of the step
[6140]368         IF( nb_trcobc > 0 ) THEN
369           if (lwp) write(numout,'(a,i5,a,i10)') '   reading OBC data for ', nb_trcobc ,' variable(s) at step ', kt
[12256]370           CALL fld_read( kt=kt, kn_fsbc=1, sd=sf_trcobc, kit=jit, pt_offset = 0.5_wp )
[6140]371         ENDIF
[9169]372         !
[6140]373         ! SURFACE boundary conditions
374         IF( nb_trcsbc > 0 ) THEN
375           if (lwp) write(numout,'(a,i5,a,i10)') '   reading SBC data for ', nb_trcsbc ,' variable(s) at step ', kt
[9169]376           CALL fld_read( kt=kt, kn_fsbc=1, sd=sf_trcsbc, kit=jit)
[6140]377         ENDIF
[9169]378         !
[6140]379         ! COASTAL boundary conditions
380         IF( nb_trccbc > 0 ) THEN
381           if (lwp) write(numout,'(a,i5,a,i10)') '   reading CBC data for ', nb_trccbc ,' variable(s) at step ', kt
[9169]382           CALL fld_read( kt=kt, kn_fsbc=1, sd=sf_trccbc, kit=jit)
[6140]383         ENDIF
[9169]384         !
[6140]385      ELSE
[9169]386         !
[12256]387         ! BDY: use pt_offset=0.5 as applied at the end of the step and fldread is referenced at the middle of the step
[6140]388         IF( nb_trcobc > 0 ) THEN
389           if (lwp) write(numout,'(a,i5,a,i10)') '   reading OBC data for ', nb_trcobc ,' variable(s) at step ', kt
[12256]390           CALL fld_read( kt=kt, kn_fsbc=1, sd=sf_trcobc, pt_offset = 0.5_wp )
[6140]391         ENDIF
[9169]392         !
[6140]393         ! SURFACE boundary conditions
394         IF( nb_trcsbc > 0 ) THEN
395           if (lwp) write(numout,'(a,i5,a,i10)') '   reading SBC data for ', nb_trcsbc ,' variable(s) at step ', kt
[9169]396           CALL fld_read( kt=kt, kn_fsbc=1, sd=sf_trcsbc )
[6140]397         ENDIF
[9169]398         !
[6140]399         ! COASTAL boundary conditions
400         IF( nb_trccbc > 0 ) THEN
401           if (lwp) write(numout,'(a,i5,a,i10)') '   reading CBC data for ', nb_trccbc ,' variable(s) at step ', kt
[9169]402           CALL fld_read( kt=kt, kn_fsbc=1, sd=sf_trccbc )
[6140]403         ENDIF
[9169]404         !
[4058]405      ENDIF
406
[7646]407      ! 2. Apply Boundary conditions data
408      !
409      DO jn = 1 , jptra
410         !
411         ! Remove river dilution for tracers with absent river load
[9169]412         IF( ln_rnf_ctl .AND. .NOT.ln_trc_cbc(jn) ) THEN
[7646]413            DO jj = 2, jpj
414               DO ji = fs_2, fs_jpim1
415                  DO jk = 1, nk_rnf(ji,jj)
416                     zrnf = (rnf(ji,jj) + rnf_b(ji,jj)) * 0.5_wp * r1_rau0 / h_rnf(ji,jj)
417                     tra(ji,jj,jk,jn) = tra(ji,jj,jk,jn)  + (trn(ji,jj,jk,jn) * zrnf)
[9169]418                  END DO
419               END DO
420            END DO
[7646]421         ENDIF
[9169]422         !
[7646]423         ! OPEN boundary conditions: trcbdy is called in trcnxt !
[9169]424         !
[7646]425         ! SURFACE boundary conditions
[9169]426         IF( ln_trc_sbc(jn) ) THEN
[7646]427            jl = n_trc_indsbc(jn)
[12110]428            sf_trcsbc(jl)%fnow(:,:,1) = MAX( rtrn, sf_trcsbc(jl)%fnow(:,:,1) ) ! avoid nedgative value due to interpolation
[7646]429            DO jj = 2, jpj
430               DO ji = fs_2, fs_jpim1   ! vector opt.
[12110]431                  zfact = 1. / ( e3t_n(ji,jj,1) * rn_sbc_time )
[7646]432                  tra(ji,jj,1,jn) = tra(ji,jj,1,jn) + rf_trsfac(jl) * sf_trcsbc(jl)%fnow(ji,jj,1) * zfact
433               END DO
434            END DO
[9169]435         ENDIF
436         !
[7646]437         ! COASTAL boundary conditions
[12110]438         IF( ( ln_rnf .OR. l_offline ) .AND. ln_trc_cbc(jn) ) THEN
439            IF( l_offline )   rn_rfact = 1._wp
[7646]440            jl = n_trc_indcbc(jn)
441            DO jj = 2, jpj
442               DO ji = fs_2, fs_jpim1   ! vector opt.
443                  DO jk = 1, nk_rnf(ji,jj)
[12110]444                     zfact = rn_rfact / ( e1e2t(ji,jj) * h_rnf(ji,jj) * rn_cbc_time ) * tmask(ji,jj,1)
[7646]445                     tra(ji,jj,jk,jn) = tra(ji,jj,jk,jn) + rf_trcfac(jl) * sf_trccbc(jl)%fnow(ji,jj,1) * zfact
[9169]446                  END DO
[7646]447               END DO
448            END DO
[9169]449         ENDIF
[7646]450         !                                                       ! ===========
451      END DO                                                     ! tracer loop
452      !                                                          ! ===========
[9124]453      IF( ln_timing )   CALL timing_stop('trc_bc')
[4058]454      !
[7646]455   END SUBROUTINE trc_bc
[4058]456
457#else
458   !!----------------------------------------------------------------------
459   !!   Dummy module                              NO 3D passive tracer data
460   !!----------------------------------------------------------------------
[12209]461   
462 CONTAINS
463
[7646]464   SUBROUTINE trc_bc_ini( ntrc )        ! Empty routine
[6140]465      INTEGER,INTENT(IN) :: ntrc                           ! number of tracers
[12327]466      WRITE(*,*) 'trc_bc_ini: You should not have seen this print! error?', ntrc
[7646]467   END SUBROUTINE trc_bc_ini
[12327]468   SUBROUTINE trc_bc( kt, jit )        ! Empty routine
469      INTEGER,INTENT(IN) :: kt, jit                           ! number of tracers
470      WRITE(*,*) 'trc_bc: You should not have seen this print! error?', kt, jit
[7646]471   END SUBROUTINE trc_bc
[4058]472#endif
473
474   !!======================================================================
475END MODULE trcbc
Note: See TracBrowser for help on using the repository browser.