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 branches/2016/dev_merge_2016/NEMOGCM/NEMO/TOP_SRC – NEMO

source: branches/2016/dev_merge_2016/NEMOGCM/NEMO/TOP_SRC/trcbc.F90 @ 7403

Last change on this file since 7403 was 7403, checked in by timgraham, 8 years ago

Merge dev_INGV_METO_merge_2016 into branch

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