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.
Changeset 10251 for branches/UKMO/dev_r5518_AMM15_package/NEMOGCM/NEMO/TOP_SRC/trcbc.F90 – NEMO

Ignore:
Timestamp:
2018-10-29T15:20:26+01:00 (5 years ago)
Author:
kingr
Message:

Rolled back to r10247 - i.e., undid merge of pkg br and 3.6_stable br

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/UKMO/dev_r5518_AMM15_package/NEMOGCM/NEMO/TOP_SRC/trcbc.F90

    r10249 r10251  
    44   !! TOP :  module for passive tracer boundary conditions 
    55   !!===================================================================== 
    6    !! History :  3.5 !  2014-04  (M. Vichi, T. Lovato)  Original 
    7    !!            3.6 !  2015-03  (T . Lovato) Revision and BDY support 
    8    !!---------------------------------------------------------------------- 
    9 #if defined key_top 
     6   !!---------------------------------------------------------------------- 
     7#if  defined key_top  
    108   !!---------------------------------------------------------------------- 
    119   !!   'key_top'                                                TOP model  
    1210   !!---------------------------------------------------------------------- 
    13    !!   trc_bc       : read and time interpolated tracer Boundary Conditions 
     11   !!   trc_dta    : read and time interpolated passive tracer data 
    1412   !!---------------------------------------------------------------------- 
    1513   USE par_trc       !  passive tracers parameters 
     
    1917   USE lib_mpp       !  MPP library 
    2018   USE fldread       !  read input fields 
    21 #if defined key_bdy 
    22    USE bdy_oce, only: nb_bdy , idx_bdy, ln_coords_file, rn_time_dmp, rn_time_dmp_out 
    23 #endif 
    2419 
    2520   IMPLICIT NONE 
     
    2924   PUBLIC   trc_bc_read    ! called in trcstp.F90 or within 
    3025 
    31    INTEGER  , SAVE, PUBLIC                             :: nb_trcobc    ! number of tracers with open BC 
    32    INTEGER  , SAVE, PUBLIC                             :: nb_trcsbc    ! number of tracers with surface BC 
    33    INTEGER  , SAVE, PUBLIC                             :: nb_trccbc    ! number of tracers with coastal BC 
     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 
    3429   INTEGER  , SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: n_trc_indobc ! index of tracer with OBC data 
    3530   INTEGER  , SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: n_trc_indsbc ! index of tracer with SBC data 
    3631   INTEGER  , SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: n_trc_indcbc ! index of tracer with CBC data 
    37    REAL(wp) , SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: rf_trsfac    ! multiplicative factor for SBC tracer values 
    38    TYPE(FLD), SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: sf_trcsbc    ! structure of data input SBC (file informations, fields read) 
    39    REAL(wp) , SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: rf_trcfac    ! multiplicative factor for CBC tracer values 
    40    TYPE(FLD), SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: sf_trccbc    ! structure of data input CBC (file informations, fields read) 
    41    REAL(wp) , SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:,:)  :: rf_trofac    ! multiplicative factor for OBCtracer values 
    42    TYPE(FLD), SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:,:), TARGET  :: sf_trcobc    ! structure of data input OBC (file informations, fields read) 
    43    TYPE(MAP_POINTER), ALLOCATABLE, DIMENSION(:,:) :: nbmap_ptr   ! array of pointers to nbmap 
     32   INTEGER  , SAVE, PUBLIC                             :: ntra_obc     ! MAX( 1, nb_trcxxx ) to avoid compilation error with bounds checking 
     33   INTEGER  , SAVE, PUBLIC                             :: ntra_sbc     ! MAX( 1, nb_trcxxx ) to avoid compilation error with bounds checking 
     34   INTEGER  , SAVE, PUBLIC                             :: ntra_cbc     ! MAX( 1, nb_trcxxx ) to avoid compilation error with bounds checking 
     35   REAL(wp) , SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: rf_trofac   ! multiplicative factor for OBCtracer values 
     36   TYPE(FLD), SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: sf_trcobc   ! structure of data input OBC (file informations, fields read) 
     37   REAL(wp) , SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: rf_trsfac   ! multiplicative factor for SBC tracer values 
     38   TYPE(FLD), SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: sf_trcsbc   ! structure of data input SBC (file informations, fields read) 
     39   REAL(wp) , SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: rf_trcfac   ! multiplicative factor for CBC tracer values 
     40   TYPE(FLD), SAVE, PUBLIC, ALLOCATABLE, DIMENSION(:)  :: sf_trccbc   ! structure of data input CBC (file informations, fields read) 
    4441 
    4542   !! * Substitutions 
    4643#  include "domzgr_substitute.h90" 
    4744   !!---------------------------------------------------------------------- 
    48    !! NEMO/OPA 3.6 , NEMO Consortium (2015) 
     45   !! NEMO/OPA 3.3 , NEMO Consortium (2010) 
    4946   !! $Id$ 
    5047   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt) 
     
    6360      ! 
    6461      INTEGER,INTENT(IN) :: ntrc                           ! number of tracers 
    65       INTEGER            :: jl, jn , ib, ibd, ii, ij, ik   ! dummy loop indices 
     62      INTEGER            :: jl, jn                         ! dummy loop indices 
    6663      INTEGER            :: ierr0, ierr1, ierr2, ierr3     ! temporary integers 
    67       INTEGER            :: ios                            ! Local integer output status for namelist read 
    68       INTEGER            :: nblen, igrd                    ! support arrays for BDY 
     64      INTEGER            ::  ios                           ! Local integer output status for namelist read 
    6965      CHARACTER(len=100) :: clndta, clntrc 
    7066      ! 
    71       CHARACTER(len=100) :: cn_dir_sbc, cn_dir_cbc, cn_dir_obc 
    72  
     67      CHARACTER(len=100) :: cn_dir 
    7368      TYPE(FLD_N), ALLOCATABLE, DIMENSION(:) :: slf_i  ! local array of namelist informations on the fields to read 
    74       TYPE(FLD_N), DIMENSION(jpmaxtrc,2) :: sn_trcobc  ! open 
    75       TYPE(FLD_N), DIMENSION(jpmaxtrc) :: sn_trcobc2   ! to read in multiple (2) open bdy 
     69      TYPE(FLD_N), DIMENSION(jpmaxtrc) :: sn_trcobc    ! open 
    7670      TYPE(FLD_N), DIMENSION(jpmaxtrc) :: sn_trcsbc    ! surface 
    7771      TYPE(FLD_N), DIMENSION(jpmaxtrc) :: sn_trccbc    ! coastal 
     
    8074      REAL(wp)   , DIMENSION(jpmaxtrc) :: rn_trcfac    ! multiplicative factor for tracer values 
    8175      !! 
    82       NAMELIST/namtrc_bc/ cn_dir_sbc, cn_dir_cbc, cn_dir_obc, sn_trcobc2, rn_trofac, sn_trcsbc, rn_trsfac, sn_trccbc, rn_trcfac 
    83 #if defined key_bdy 
    84       NAMELIST/namtrc_bdy/ cn_trc_dflt, cn_trc, nn_trcdmp_bdy 
    85 #endif 
     76      NAMELIST/namtrc_bc/ cn_dir, sn_trcobc, rn_trofac, sn_trcsbc, rn_trsfac, sn_trccbc, rn_trcfac  
    8677      !!---------------------------------------------------------------------- 
    8778      IF( nn_timing == 1 )  CALL timing_start('trc_bc_init') 
    8879      ! 
    89       IF( lwp ) THEN 
    90          WRITE(numout,*) ' ' 
    91          WRITE(numout,*) 'trc_bc_init : Tracers Boundary Conditions (BC)' 
    92          WRITE(numout,*) '~~~~~~~~~~~ ' 
    93       ENDIF 
    9480      !  Initialisation and local array allocation 
    9581      ierr0 = 0  ;  ierr1 = 0  ;  ierr2 = 0  ;  ierr3 = 0   
     
    121107      n_trc_indcbc(:) = 0 
    122108      ! 
    123       ! Read Boundary Conditions Namelists 
     109      DO jn = 1, ntrc 
     110         IF( ln_trc_obc(jn) ) THEN 
     111             nb_trcobc       = nb_trcobc + 1  
     112             n_trc_indobc(jn) = nb_trcobc  
     113         ENDIF 
     114         IF( ln_trc_sbc(jn) ) THEN 
     115             nb_trcsbc       = nb_trcsbc + 1 
     116             n_trc_indsbc(jn) = nb_trcsbc 
     117         ENDIF 
     118         IF( ln_trc_cbc(jn) ) THEN 
     119             nb_trccbc       = nb_trccbc + 1 
     120             n_trc_indcbc(jn) = nb_trccbc 
     121         ENDIF 
     122      ENDDO 
     123      ntra_obc = MAX( 1, nb_trcobc )   ! To avoid compilation error with bounds checking 
     124      IF( lwp ) WRITE(numout,*) ' ' 
     125      IF( lwp ) WRITE(numout,*) ' Number of passive tracers to be initialized with open boundary data :', nb_trcobc 
     126      IF( lwp ) WRITE(numout,*) ' ' 
     127      ntra_sbc = MAX( 1, nb_trcsbc )   ! To avoid compilation error with bounds checking 
     128      IF( lwp ) WRITE(numout,*) ' ' 
     129      IF( lwp ) WRITE(numout,*) ' Number of passive tracers to be initialized with surface boundary data :', nb_trcsbc 
     130      IF( lwp ) WRITE(numout,*) ' ' 
     131      ntra_cbc = MAX( 1, nb_trccbc )   ! To avoid compilation error with bounds checking 
     132      IF( lwp ) WRITE(numout,*) ' ' 
     133      IF( lwp ) WRITE(numout,*) ' Number of passive tracers to be initialized with coastal boundary data :', nb_trccbc 
     134      IF( lwp ) WRITE(numout,*) ' ' 
     135 
    124136      REWIND( numnat_ref )              ! Namelist namtrc_bc in reference namelist : Passive tracer data structure 
    125137      READ  ( numnat_ref, namtrc_bc, IOSTAT = ios, ERR = 901) 
     
    127139 
    128140      REWIND( numnat_cfg )              ! Namelist namtrc_bc in configuration namelist : Passive tracer data structure 
    129 #if defined key_bdy 
    130       DO ib = 1, nb_bdy 
    131 #endif 
    132         READ  ( numnat_cfg, namtrc_bc, IOSTAT = ios, ERR = 902 ) 
    133 902     IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtrc_bc in configuration namelist', lwp ) 
    134         IF(lwm) WRITE ( numont, namtrc_bc ) 
    135 #if defined key_bdy 
    136         sn_trcobc(:,ib)=sn_trcobc2(:) 
    137       ENDDO 
    138 #endif 
    139  
    140 #if defined key_bdy 
    141       REWIND( numnat_ref )              ! Namelist namtrc_bc in reference namelist : Passive tracer data structure 
    142       READ  ( numnat_ref, namtrc_bdy, IOSTAT = ios, ERR = 903) 
    143 903   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtrc_bdy in reference namelist', lwp ) 
    144  
    145       REWIND( numnat_cfg )              ! Namelist namtrc_bc in configuration namelist : Passive tracer data structure 
    146       READ  ( numnat_cfg, namtrc_bdy, IOSTAT = ios, ERR = 904 ) 
    147 904   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtrc_bdy in configuration namelist', lwp ) 
    148       IF(lwm) WRITE ( numont, namtrc_bdy ) 
    149       ! setup up preliminary informations for BDY structure 
    150       DO jn = 1, ntrc 
    151          DO ib = 1, nb_bdy 
    152             ! Set type of obc in BDY data structure (around here we may plug user override of obc type from nml) 
    153             IF ( ln_trc_obc(jn) ) THEN 
    154                trcdta_bdy(jn,ib)%cn_obc = TRIM( cn_trc(ib) ) 
    155             ELSE 
    156                trcdta_bdy(jn,ib)%cn_obc = TRIM( cn_trc_dflt(ib) ) 
    157             ENDIF 
    158             ! set damping use in BDY data structure 
    159             trcdta_bdy(jn,ib)%dmp = .false. 
    160             IF(nn_trcdmp_bdy(ib) .EQ. 1 .AND. ln_trc_obc(jn) ) trcdta_bdy(jn,ib)%dmp = .true. 
    161             IF(nn_trcdmp_bdy(ib) .EQ. 2 ) trcdta_bdy(jn,ib)%dmp = .true. 
    162             IF(trcdta_bdy(jn,ib)%cn_obc == 'frs' .AND. nn_trcdmp_bdy(ib) .NE. 0 )  & 
    163                 & CALL ctl_stop( 'Use FRS OR relaxation' ) 
    164             IF (nn_trcdmp_bdy(ib) .LT. 0 .OR. nn_trcdmp_bdy(ib) .GT. 2)            & 
    165                 & CALL ctl_stop( 'Not a valid option for nn_trcdmp_bdy. Allowed: 0,1,2.' ) 
     141      READ  ( numnat_cfg, namtrc_bc, IOSTAT = ios, ERR = 902 ) 
     142902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtrc_bc in configuration namelist', lwp ) 
     143      IF(lwm) WRITE ( numont, namtrc_bc ) 
     144 
     145      ! print some information for each  
     146      IF( lwp ) THEN 
     147         DO jn = 1, ntrc 
     148            IF( ln_trc_obc(jn) )  THEN     
     149               clndta = TRIM( sn_trcobc(jn)%clvar )  
     150               IF(lwp) WRITE(numout,*) 'Preparing to read OBC data file for passive tracer number :', jn, ' name : ', clndta, &  
     151               &               ' multiplicative factor : ', rn_trofac(jn) 
     152            ENDIF 
     153            IF( ln_trc_sbc(jn) )  THEN     
     154               clndta = TRIM( sn_trcsbc(jn)%clvar )  
     155               IF(lwp) WRITE(numout,*) 'Preparing to read SBC data file for passive tracer number :', jn, ' name : ', clndta, &  
     156               &               ' multiplicative factor : ', rn_trsfac(jn) 
     157            ENDIF 
     158            IF( ln_trc_cbc(jn) )  THEN     
     159               clndta = TRIM( sn_trccbc(jn)%clvar )  
     160               IF(lwp) WRITE(numout,*) 'Preparing to read CBC data file for passive tracer number :', jn, ' name : ', clndta, &  
     161               &               ' multiplicative factor : ', rn_trcfac(jn) 
     162            ENDIF 
     163         END DO 
     164      ENDIF 
     165      ! 
     166      ! The following code is written this way to reduce memory usage and repeated for each boundary data 
     167      ! MAV: note that this is just a placeholder and the dimensions must be changed according to  
     168      !      what will be done with BDY. A new structure will probably need to be included 
     169      ! 
     170      ! OPEN Lateral boundary conditions 
     171      IF( nb_trcobc > 0 ) THEN       !  allocate only if the number of tracer to initialise is greater than zero 
     172         ALLOCATE( sf_trcobc(nb_trcobc), rf_trofac(nb_trcobc), STAT=ierr1 ) 
     173         IF( ierr1 > 0 ) THEN 
     174            CALL ctl_stop( 'trc_bc_init: unable to allocate  sf_trcobc structure' )   ;   RETURN 
     175         ENDIF 
     176         ! 
     177         DO jn = 1, ntrc 
     178            IF( ln_trc_obc(jn) ) THEN      ! update passive tracers arrays with input data read from file 
     179               jl = n_trc_indobc(jn) 
     180               slf_i(jl)    = sn_trcobc(jn) 
     181               rf_trofac(jl) = rn_trofac(jn) 
     182                                            ALLOCATE( sf_trcobc(jl)%fnow(jpi,jpj,jpk)   , STAT=ierr2 ) 
     183               IF( sn_trcobc(jn)%ln_tint )  ALLOCATE( sf_trcobc(jl)%fdta(jpi,jpj,jpk,2) , STAT=ierr3 ) 
     184               IF( ierr2 + ierr3 > 0 ) THEN 
     185                 CALL ctl_stop( 'trc_bc_init : unable to allocate passive tracer OBC data arrays' )   ;   RETURN 
     186               ENDIF 
     187            ENDIF 
     188            !    
    166189         ENDDO 
    167       ENDDO 
    168  
    169 #else 
    170       ! Force all tracers OBC to false if bdy not used 
    171       ln_trc_obc = .false. 
    172 #endif 
    173       ! compose BC data indexes 
    174       DO jn = 1, ntrc 
    175          IF( ln_trc_obc(jn) ) THEN 
    176              nb_trcobc       = nb_trcobc + 1  ; n_trc_indobc(jn) = nb_trcobc 
    177          ENDIF 
    178          IF( ln_trc_sbc(jn) ) THEN 
    179              nb_trcsbc       = nb_trcsbc + 1  ; n_trc_indsbc(jn) = nb_trcsbc 
    180          ENDIF 
    181          IF( ln_trc_cbc(jn) ) THEN 
    182              nb_trccbc       = nb_trccbc + 1  ; n_trc_indcbc(jn) = nb_trccbc 
    183          ENDIF 
    184       ENDDO 
    185  
    186       ! Print summmary of Boundary Conditions 
    187       IF( lwp ) THEN 
    188          WRITE(numout,*) ' ' 
    189          WRITE(numout,'(a,i3)') '   Total tracers to be initialized with SURFACE BCs data:', nb_trcsbc 
    190          IF ( nb_trcsbc > 0 ) THEN 
    191             WRITE(numout,*) '   #trc        NAME        Boundary     Mult.Fact. ' 
    192             DO jn = 1, ntrc 
    193                IF ( ln_trc_sbc(jn) ) WRITE(numout,9001) jn, TRIM( sn_trcsbc(jn)%clvar ), 'SBC', rn_trsfac(jn) 
    194             ENDDO 
    195          ENDIF 
    196          WRITE(numout,'(2a)') '   SURFACE BC data repository : ', TRIM(cn_dir_sbc) 
    197  
    198          WRITE(numout,*) ' ' 
    199          WRITE(numout,'(a,i3)') '   Total tracers to be initialized with COASTAL BCs data:', nb_trccbc 
    200          IF ( nb_trccbc > 0 ) THEN 
    201             WRITE(numout,*) '   #trc        NAME        Boundary     Mult.Fact. ' 
    202             DO jn = 1, ntrc 
    203                IF ( ln_trc_cbc(jn) ) WRITE(numout, 9001) jn, TRIM( sn_trccbc(jn)%clvar ), 'CBC', rn_trcfac(jn) 
    204             ENDDO 
    205          ENDIF 
    206          WRITE(numout,'(2a)') '   COASTAL BC data repository : ', TRIM(cn_dir_cbc) 
    207  
    208          WRITE(numout,*) ' ' 
    209          WRITE(numout,'(a,i3)') '   Total tracers to be initialized with OPEN BCs data:', nb_trcobc 
    210 #if defined key_bdy 
    211          IF ( nb_trcobc > 0 ) THEN 
    212             WRITE(numout,*) '   #trc        NAME        Boundary     Mult.Fact.   OBC Settings' 
    213             DO jn = 1, ntrc 
    214                DO ib = 1, nb_bdy 
    215                  IF ( ln_trc_obc(jn) )  WRITE(numout, 9001) jn, TRIM( sn_trcobc(jn,ib)%clvar ), 'OBC', rn_trofac(jn), (trcdta_bdy(jn,ib)%cn_obc) 
    216                ENDDO 
    217                !IF ( ln_trc_obc(jn) )  WRITE(numout, 9001) jn, TRIM( sn_trcobc(jn,ib)%clvar ), 'OBC', rn_trofac(jn), (trcdta_bdy(jn,ib)%cn_obc,ib=1,nb_bdy) 
    218                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) 
    219             ENDDO 
    220             WRITE(numout,*) ' ' 
    221             DO ib = 1, nb_bdy 
    222                 IF (nn_trcdmp_bdy(ib) .EQ. 0) WRITE(numout,9003) '   Boundary ',ib,' -> NO damping of tracers' 
    223                 IF (nn_trcdmp_bdy(ib) .EQ. 1) WRITE(numout,9003) '   Boundary ',ib,' -> damping ONLY for tracers with external data provided' 
    224                 IF (nn_trcdmp_bdy(ib) .EQ. 2) WRITE(numout,9003) '   Boundary ',ib,' -> damping of ALL tracers' 
    225                 IF (nn_trcdmp_bdy(ib) .GT. 0) THEN 
    226                    WRITE(numout,9003) '     USE damping parameters from nambdy for boundary ', ib,' : ' 
    227                    WRITE(numout,'(a,f10.2,a)') '     - Inflow damping time scale  : ',rn_time_dmp(ib),' days' 
    228                    WRITE(numout,'(a,f10.2,a)') '     - Outflow damping time scale : ',rn_time_dmp_out(ib),' days' 
    229                 ENDIF 
    230             ENDDO 
    231          ENDIF 
    232 #endif 
    233          WRITE(numout,'(2a)') '   OPEN BC data repository : ', TRIM(cn_dir_obc) 
    234       ENDIF 
    235 9001  FORMAT(2x,i5, 3x, a15, 3x, a5, 6x, e11.3, 4x, 10a13) 
    236 9002  FORMAT(2x,i5, 3x, a41, 3x, 10a13) 
    237 9003  FORMAT(a, i5, a) 
    238  
    239       ! 
    240 #if defined key_bdy 
    241       ! OPEN Lateral boundary conditions 
    242       IF( nb_trcobc > 0 ) THEN  
    243          ALLOCATE ( sf_trcobc(nb_trcobc,nb_bdy), rf_trofac(nb_trcobc,nb_bdy), nbmap_ptr(nb_trcobc,nb_bdy), STAT=ierr1 ) 
    244          IF( ierr1 > 0 ) THEN 
    245             CALL ctl_stop( 'trc_bc_init: unable to allocate sf_trcobc structure' )   ;   RETURN 
    246          ENDIF 
    247  
    248          igrd = 1                       ! Everything is at T-points here 
    249  
    250          DO ib = 1, nb_bdy 
    251             DO jn = 1, ntrc 
    252  
    253                nblen = idx_bdy(ib)%nblen(igrd) 
    254  
    255                IF ( ln_trc_obc(jn) ) THEN 
    256                ! Initialise from external data 
    257                   jl = n_trc_indobc(jn) 
    258                   slf_i(jl)    = sn_trcobc(jn,ib) 
    259                   rf_trofac(jl,ib) = rn_trofac(jn) 
    260                                                ALLOCATE( sf_trcobc(jl,ib)%fnow(nblen,1,jpk)   , STAT=ierr2 ) 
    261                   IF( sn_trcobc(jn,ib)%ln_tint )  ALLOCATE( sf_trcobc(jl,ib)%fdta(nblen,1,jpk,2) , STAT=ierr3 ) 
    262                   IF( ierr2 + ierr3 > 0 ) THEN 
    263                     CALL ctl_stop( 'trc_bc_init : unable to allocate passive tracer OBC data arrays' )   ;   RETURN 
    264                   ENDIF 
    265                   trcdta_bdy(jn,ib)%trc => sf_trcobc(jl,ib)%fnow(:,1,:) 
    266                   trcdta_bdy(jn,ib)%rn_fac = rf_trofac(jl,ib) 
    267                   ! create OBC mapping array 
    268                   nbmap_ptr(jl,ib)%ptr => idx_bdy(ib)%nbmap(:,igrd) 
    269                   nbmap_ptr(jl,ib)%ll_unstruc = ln_coords_file(igrd) 
    270                ELSE 
    271                ! Initialise obc arrays from initial conditions 
    272                   ALLOCATE ( trcdta_bdy(jn,ib)%trc(nblen,jpk) ) 
    273                   DO ibd = 1, nblen 
    274                      DO ik = 1, jpkm1 
    275                         ii = idx_bdy(ib)%nbi(ibd,igrd) 
    276                         ij = idx_bdy(ib)%nbj(ibd,igrd) 
    277                         trcdta_bdy(jn,ib)%trc(ibd,ik) = trn(ii,ij,ik,jn) * tmask(ii,ij,ik) 
    278                      END DO 
    279                   END DO 
    280                   trcdta_bdy(jn,ib)%rn_fac = 1._wp 
    281                ENDIF 
    282             ENDDO 
    283             CALL fld_fill( sf_trcobc(:,ib), slf_i, cn_dir_obc, 'trc_bc_init', 'Passive tracer OBC data', 'namtrc_bc' ) 
    284          ENDDO 
    285  
    286       ENDIF 
    287 #endif 
     190         !                         ! fill sf_trcdta with slf_i and control print 
     191         CALL fld_fill( sf_trcobc, slf_i, cn_dir, 'trc_bc_init', 'Passive tracer OBC data', 'namtrc_bc' ) 
     192         ! 
     193      ENDIF 
     194      ! 
    288195      ! SURFACE Boundary conditions 
    289196      IF( nb_trcsbc > 0 ) THEN       !  allocate only if the number of tracer to initialise is greater than zero 
     
    307214         ENDDO 
    308215         !                         ! fill sf_trcsbc with slf_i and control print 
    309          CALL fld_fill( sf_trcsbc, slf_i, cn_dir_sbc, 'trc_bc_init', 'Passive tracer SBC data', 'namtrc_bc' ) 
     216         CALL fld_fill( sf_trcsbc, slf_i, cn_dir, 'trc_bc_init', 'Passive tracer SBC data', 'namtrc_bc' ) 
    310217         ! 
    311218      ENDIF 
     
    332239         ENDDO 
    333240         !                         ! fill sf_trccbc with slf_i and control print 
    334          CALL fld_fill( sf_trccbc, slf_i, cn_dir_cbc, 'trc_bc_init', 'Passive tracer CBC data', 'namtrc_bc' ) 
     241         CALL fld_fill( sf_trccbc, slf_i, cn_dir, 'trc_bc_init', 'Passive tracer CBC data', 'namtrc_bc' ) 
    335242         ! 
    336243      ENDIF 
     
    342249 
    343250 
    344    SUBROUTINE trc_bc_read(kt, jit) 
     251   SUBROUTINE trc_bc_read(kt) 
    345252      !!---------------------------------------------------------------------- 
    346253      !!                   ***  ROUTINE trc_bc_init  *** 
     
    357264      !! * Arguments 
    358265      INTEGER, INTENT( in ) ::   kt      ! ocean time-step index 
    359       INTEGER, INTENT( in ), OPTIONAL ::   jit   ! subcycle time-step index (for timesplitting option) 
    360       INTEGER :: ib 
     266 
    361267      !!--------------------------------------------------------------------- 
    362268      ! 
    363269      IF( nn_timing == 1 )  CALL timing_start('trc_bc_read') 
    364270 
    365       IF( kt == nit000 .AND. lwp) THEN 
    366          WRITE(numout,*) 
    367          WRITE(numout,*) 'trc_bc_read : Surface boundary conditions for passive tracers.' 
    368          WRITE(numout,*) '~~~~~~~~~~~ ' 
    369       ENDIF 
    370  
    371       IF ( PRESENT(jit) ) THEN  
    372  
    373 #ifdef key_bdy 
    374          ! OPEN boundary conditions (use time_offset=+1 as they are applied at the end of the step) 
    375          IF( nb_trcobc > 0 ) THEN 
    376            if (lwp) write(numout,'(a,i5,a,i10)') '   reading OBC data for ', nb_trcobc ,' variable(s) at step ', kt 
    377            DO ib = 1,nb_bdy 
    378              CALL fld_read(kt=kt, kn_fsbc=1, sd=sf_trcobc(:,ib), map=nbmap_ptr(:,ib), kit=jit, kt_offset=+1) 
    379            ENDDO 
    380          ENDIF 
    381 #endif 
    382  
    383          ! SURFACE boundary conditions 
    384          IF( nb_trcsbc > 0 ) THEN 
    385            if (lwp) write(numout,'(a,i5,a,i10)') '   reading SBC data for ', nb_trcsbc ,' variable(s) at step ', kt 
    386            CALL fld_read(kt=kt, kn_fsbc=1, sd=sf_trcsbc, kit=jit) 
    387          ENDIF 
    388  
    389          ! COASTAL boundary conditions 
    390          IF( nb_trccbc > 0 ) THEN 
    391            if (lwp) write(numout,'(a,i5,a,i10)') '   reading CBC data for ', nb_trccbc ,' variable(s) at step ', kt 
    392            CALL fld_read(kt=kt, kn_fsbc=1, sd=sf_trccbc, kit=jit) 
    393          ENDIF 
    394  
    395       ELSE 
    396  
    397 #ifdef key_bdy 
    398          ! OPEN boundary conditions (use time_offset=+1 as they are applied at the end of the step) 
    399          IF( nb_trcobc > 0 ) THEN 
    400            if (lwp) write(numout,'(a,i5,a,i10)') '   reading OBC data for ', nb_trcobc ,' variable(s) at step ', kt 
    401            DO ib = 1,nb_bdy 
    402              CALL fld_read(kt=kt, kn_fsbc=1, sd=sf_trcobc(:,ib), map=nbmap_ptr(:,ib), kt_offset=+1) 
    403            ENDDO 
    404          ENDIF 
    405 #endif 
    406  
    407          ! SURFACE boundary conditions 
    408          IF( nb_trcsbc > 0 ) THEN 
    409            if (lwp) write(numout,'(a,i5,a,i10)') '   reading SBC data for ', nb_trcsbc ,' variable(s) at step ', kt 
    410            CALL fld_read(kt=kt, kn_fsbc=1, sd=sf_trcsbc) 
    411          ENDIF 
    412  
    413          ! COASTAL boundary conditions 
    414          IF( nb_trccbc > 0 ) THEN 
    415            if (lwp) write(numout,'(a,i5,a,i10)') '   reading CBC data for ', nb_trccbc ,' variable(s) at step ', kt 
    416            CALL fld_read(kt=kt, kn_fsbc=1, sd=sf_trccbc) 
    417          ENDIF 
    418  
    419       ENDIF 
    420  
     271      IF( kt == nit000 ) THEN 
     272         IF(lwp) WRITE(numout,*) 
     273         IF(lwp) WRITE(numout,*) 'trc_bc_read : Surface boundary conditions for passive tracers.' 
     274         IF(lwp) WRITE(numout,*) '~~~~~~~ ' 
     275      ENDIF 
     276 
     277      ! OPEN boundary conditions: DOES NOT WORK. Waiting for stable BDY 
     278      IF( nb_trcobc > 0 ) THEN 
     279        if (lwp) write(numout,'(a,i5,a,i5)') '   reading OBC data for ', nb_trcobc ,' variables at step ', kt 
     280        CALL fld_read(kt,1,sf_trcobc) 
     281        ! vertical interpolation on s-grid and partial step to be added 
     282      ENDIF 
     283 
     284      ! SURFACE boundary conditions        
     285      IF( nb_trcsbc > 0 ) THEN 
     286        if (lwp) write(numout,'(a,i5,a,i5)') '   reading SBC data for ', nb_trcsbc ,' variables at step ', kt 
     287        CALL fld_read(kt,1,sf_trcsbc) 
     288      ENDIF 
     289 
     290      ! COASTAL boundary conditions        
     291      IF( nb_trccbc > 0 ) THEN 
     292        if (lwp) write(numout,'(a,i5,a,i5)') '   reading CBC data for ', nb_trccbc ,' variables at step ', kt 
     293        CALL fld_read(kt,1,sf_trccbc) 
     294      ENDIF    
    421295      ! 
    422296      IF( nn_timing == 1 )  CALL timing_stop('trc_bc_read') 
     
    429303   !!---------------------------------------------------------------------- 
    430304CONTAINS 
    431  
    432    SUBROUTINE trc_bc_init( ntrc )        ! Empty routine 
    433       INTEGER,INTENT(IN) :: ntrc                           ! number of tracers 
    434       WRITE(*,*) 'trc_bc_init: You should not have seen this print! error?', kt 
    435    END SUBROUTINE trc_bc_init 
    436  
    437305   SUBROUTINE trc_bc_read( kt )        ! Empty routine 
    438306      WRITE(*,*) 'trc_bc_read: You should not have seen this print! error?', kt 
Note: See TracChangeset for help on using the changeset viewer.