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/2020/dev_r12563_ASINTER-06_ABL_improvement/src/TOP – NEMO

source: NEMO/branches/2020/dev_r12563_ASINTER-06_ABL_improvement/src/TOP/trcbc.F90 @ 12808

Last change on this file since 12808 was 12489, checked in by davestorkey, 4 years ago

Preparation for new timestepping scheme #2390.
Main changes:

  1. Initial euler timestep now handled in stp and not in TRA/DYN routines.
  2. Renaming of all timestep parameters. In summary, the namelist parameter is now rn_Dt and the current timestep is rDt (and rDt_ice, rDt_trc etc).
  3. Renaming of a few miscellaneous parameters, eg. atfp -> rn_atfp (namelist parameter used everywhere) and rau0 -> rho0.

This version gives bit-comparable results to the previous version of the trunk.

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