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 12205 for NEMO/branches/2019/dev_r11943_MERGE_2019/src – NEMO

Ignore:
Timestamp:
2019-12-12T11:52:50+01:00 (4 years ago)
Author:
acc
Message:

2019/dev_r11943_MERGE_2019: Merge in dev_r11879_ENHANCE-05_SimonM-Harmonic_Analysis. SETTE tested

Location:
NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE
Files:
6 deleted
10 edited
3 copied

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/BDY/bdydta.F90

    r12145 r12205  
    2323   USE phycst         ! physical constants 
    2424   USE sbcapr         ! atmospheric pressure forcing 
    25    USE sbctide        ! Tidal forcing or not 
     25   USE tide_mod, ONLY: ln_tide ! tidal forcing 
    2626   USE bdy_oce        ! ocean open boundary conditions   
    2727   USE bdytides       ! tidal forcing at boundaries 
     
    7575CONTAINS 
    7676 
    77    SUBROUTINE bdy_dta( kt, Kmm, kit, kt_offset ) 
     77   SUBROUTINE bdy_dta( kt, Kmm, kt_offset ) 
    7878      !!---------------------------------------------------------------------- 
    7979      !!                   ***  SUBROUTINE bdy_dta  *** 
     
    8686      INTEGER, INTENT(in)           ::   kt           ! ocean time-step index  
    8787      INTEGER, INTENT(in)           ::   Kmm          ! ocean time level index 
    88       INTEGER, INTENT(in), OPTIONAL ::   kit          ! subcycle time-step index (for timesplitting option) 
    89       INTEGER, INTENT(in), OPTIONAL ::   kt_offset    ! time offset in units of timesteps. NB. if kit 
     88      INTEGER, INTENT(in), OPTIONAL ::   kt_offset    ! time offset in units of timesteps 
    9089      !                                               ! is present then units = subcycle timesteps. 
    9190      !                                               ! kt_offset = 0 => get data at "now" time level 
     
    106105      ! Initialise data arrays once for all from initial conditions where required 
    107106      !--------------------------------------------------------------------------- 
    108       IF( kt == nit000 .AND. .NOT.PRESENT(kit) ) THEN 
     107      IF( kt == nit000 ) THEN 
    109108 
    110109         ! Calculate depth-mean currents 
     
    217216         ! read/update all bdy data 
    218217         ! ------------------------ 
    219          CALL fld_read( kt, 1, bf_alias, kit = kit, kt_offset = kt_offset ) 
     218         CALL fld_read( kt, 1, bf_alias, kt_offset = kt_offset ) 
    220219         ! apply some corrections in some specific cases... 
    221220         ! -------------------------------------------------- 
     
    275274            END DO 
    276275         ENDIF   ! ltotvel 
    277  
    278          ! update tidal harmonic forcing 
    279          IF( PRESENT(kit) .AND. nn_dyn2d_dta(jbdy) .GE. 2 ) THEN 
    280             CALL bdytide_update( kt = kt, idx = idx_bdy(jbdy), dta = dta_alias, td = tides(jbdy),   &  
    281                &                 kit = kit, kt_offset = kt_offset ) 
    282          ENDIF 
    283276 
    284277         !  atm surface pressure : add inverted barometer effect to ssh if it was read 
  • NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/BDY/bdyini.F90

    r11960 r12205  
    2222   USE bdydta         ! open boundary cond. setting   (bdy_dta_init routine) 
    2323   USE bdytides       ! open boundary cond. setting   (bdytide_init routine) 
    24    USE sbctide        ! Tidal forcing or not 
     24   USE tide_mod, ONLY: ln_tide ! tidal forcing 
    2525   USE phycst   , ONLY: rday 
    2626   ! 
  • NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/BDY/bdytides.F90

    r11960 r12205  
    1818   USE phycst         ! physical constants 
    1919   USE bdy_oce        ! ocean open boundary conditions 
    20    USE tideini        !  
     20   USE tide_mod       !  
    2121   USE daymod         ! calendar 
    2222   ! 
     
    3030 
    3131   PUBLIC   bdytide_init     ! routine called in bdy_init 
    32    PUBLIC   bdytide_update   ! routine called in bdy_dta 
    3332   PUBLIC   bdy_dta_tides    ! routine called in dyn_spg_ts 
    3433 
     
    4544   TYPE(OBC_DATA)  , PUBLIC, DIMENSION(jp_bdy) :: dta_bdy_s  !: bdy external data (slow component) 
    4645 
     46   INTEGER ::   kt_tide 
     47 
    4748   !!---------------------------------------------------------------------- 
    4849   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
     
    6465      CHARACTER(len=80)                         ::   filtide             !: Filename root for tidal input files 
    6566      LOGICAL                                   ::   ln_bdytide_2ddta    !: If true, read 2d harmonic data 
    66       LOGICAL                                   ::   ln_bdytide_conj     !: If true, assume complex conjugate tidal data 
    6767      !! 
    6868      INTEGER                                   ::   ib_bdy, itide, ib   !: dummy loop indices 
     
    7979      TYPE(TIDES_DATA),  POINTER                ::   td                  !: local short cut    
    8080      !! 
    81       NAMELIST/nambdy_tide/filtide, ln_bdytide_2ddta, ln_bdytide_conj 
     81      NAMELIST/nambdy_tide/filtide, ln_bdytide_2ddta 
    8282      !!---------------------------------------------------------------------- 
    8383      ! 
     
    119119            IF(lwp) WRITE(numout,*) '          Namelist nambdy_tide : tidal harmonic forcing at open boundaries' 
    120120            IF(lwp) WRITE(numout,*) '             read tidal data in 2d files: ', ln_bdytide_2ddta 
    121             IF(lwp) WRITE(numout,*) '             assume complex conjugate   : ', ln_bdytide_conj 
    122121            IF(lwp) WRITE(numout,*) '             Number of tidal components to read: ', nb_harmo 
    123122            IF(lwp) THEN  
    124123                    WRITE(numout,*) '             Tidal components: '  
    125124               DO itide = 1, nb_harmo 
    126                   WRITE(numout,*)  '                 ', Wave(ntide(itide))%cname_tide  
     125                  WRITE(numout,*)  '                 ', tide_harmonics(itide)%cname_tide  
    127126               END DO 
    128127            ENDIF  
     
    165164               igrd = 1                       ! Everything is at T-points here 
    166165               DO itide = 1, nb_harmo 
    167                   CALL iom_get( inum, jpdom_autoglo, TRIM(Wave(ntide(itide))%cname_tide)//'_z1', ztr(:,:) ) 
    168                   CALL iom_get( inum, jpdom_autoglo, TRIM(Wave(ntide(itide))%cname_tide)//'_z2', zti(:,:) )  
     166                  CALL iom_get( inum, jpdom_autoglo, TRIM(tide_harmonics(itide)%cname_tide)//'_z1', ztr(:,:) ) 
     167                  CALL iom_get( inum, jpdom_autoglo, TRIM(tide_harmonics(itide)%cname_tide)//'_z2', zti(:,:) )  
    169168                  DO ib = 1, ilen0(igrd) 
    170169                     ii = idx_bdy(ib_bdy)%nbi(ib,igrd) 
     
    182181               igrd = 2                       ! Everything is at U-points here 
    183182               DO itide = 1, nb_harmo 
    184                   CALL iom_get  ( inum, jpdom_autoglo, TRIM(Wave(ntide(itide))%cname_tide)//'_u1', ztr(:,:) ) 
    185                   CALL iom_get  ( inum, jpdom_autoglo, TRIM(Wave(ntide(itide))%cname_tide)//'_u2', zti(:,:) ) 
     183                  CALL iom_get  ( inum, jpdom_autoglo, TRIM(tide_harmonics(itide)%cname_tide)//'_u1', ztr(:,:) ) 
     184                  CALL iom_get  ( inum, jpdom_autoglo, TRIM(tide_harmonics(itide)%cname_tide)//'_u2', zti(:,:) ) 
    186185                  DO ib = 1, ilen0(igrd) 
    187186                     ii = idx_bdy(ib_bdy)%nbi(ib,igrd) 
     
    199198               igrd = 3                       ! Everything is at V-points here 
    200199               DO itide = 1, nb_harmo 
    201                   CALL iom_get  ( inum, jpdom_autoglo, TRIM(Wave(ntide(itide))%cname_tide)//'_v1', ztr(:,:) ) 
    202                   CALL iom_get  ( inum, jpdom_autoglo, TRIM(Wave(ntide(itide))%cname_tide)//'_v2', zti(:,:) ) 
     200                  CALL iom_get  ( inum, jpdom_autoglo, TRIM(tide_harmonics(itide)%cname_tide)//'_v1', ztr(:,:) ) 
     201                  CALL iom_get  ( inum, jpdom_autoglo, TRIM(tide_harmonics(itide)%cname_tide)//'_v2', zti(:,:) ) 
    203202                  DO ib = 1, ilen0(igrd) 
    204203                     ii = idx_bdy(ib_bdy)%nbi(ib,igrd) 
     
    224223               DO itide = 1, nb_harmo 
    225224                  !                                                              ! SSH fields 
    226                   clfile = TRIM(filtide)//TRIM(Wave(ntide(itide))%cname_tide)//'_grid_T.nc' 
     225                  clfile = TRIM(filtide)//TRIM(tide_harmonics(itide)%cname_tide)//'_grid_T.nc' 
    227226                  CALL iom_open( clfile, inum ) 
    228227                  CALL fld_map( inum, 'z1' , dta_read(1:ilen0(1),1:1,1:1) , 1, idx_bdy(ib_bdy)%nbmap(:,1) ) 
     
    232231                  CALL iom_close( inum ) 
    233232                  !                                                              ! U fields 
    234                   clfile = TRIM(filtide)//TRIM(Wave(ntide(itide))%cname_tide)//'_grid_U.nc' 
     233                  clfile = TRIM(filtide)//TRIM(tide_harmonics(itide)%cname_tide)//'_grid_U.nc' 
    235234                  CALL iom_open( clfile, inum ) 
    236235                  CALL fld_map( inum, 'u1' , dta_read(1:ilen0(2),1:1,1:1) , 1, idx_bdy(ib_bdy)%nbmap(:,2) ) 
     
    240239                  CALL iom_close( inum ) 
    241240                  !                                                              ! V fields 
    242                   clfile = TRIM(filtide)//TRIM(Wave(ntide(itide))%cname_tide)//'_grid_V.nc' 
     241                  clfile = TRIM(filtide)//TRIM(tide_harmonics(itide)%cname_tide)//'_grid_V.nc' 
    243242                  CALL iom_open( clfile, inum ) 
    244243                  CALL fld_map( inum, 'v1' , dta_read(1:ilen0(3),1:1,1:1) , 1, idx_bdy(ib_bdy)%nbmap(:,3) ) 
     
    254253            ENDIF ! ln_bdytide_2ddta=.true. 
    255254            ! 
    256             IF( ln_bdytide_conj ) THEN    ! assume complex conjugate in data files 
    257                td%ssh0(:,:,2) = - td%ssh0(:,:,2) 
    258                td%u0  (:,:,2) = - td%u0  (:,:,2) 
    259                td%v0  (:,:,2) = - td%v0  (:,:,2) 
    260             ENDIF 
    261             ! 
    262255            ! Allocate slow varying data in the case of time splitting: 
    263256            ! Do it anyway because at this stage knowledge of free surface scheme is unknown 
     
    274267      ! 
    275268   END SUBROUTINE bdytide_init 
    276  
    277  
    278    SUBROUTINE bdytide_update( kt, idx, dta, td, kit, kt_offset ) 
    279       !!---------------------------------------------------------------------- 
    280       !!                 ***  SUBROUTINE bdytide_update  *** 
    281       !!                 
    282       !! ** Purpose : - Add tidal forcing to ssh, u2d and v2d OBC data arrays.  
    283       !!                 
    284       !!---------------------------------------------------------------------- 
    285       INTEGER          , INTENT(in   ) ::   kt          ! Main timestep counter 
    286       TYPE(OBC_INDEX)  , INTENT(in   ) ::   idx         ! OBC indices 
    287       TYPE(OBC_DATA)   , INTENT(inout) ::   dta         ! OBC external data 
    288       TYPE(TIDES_DATA) , INTENT(inout) ::   td          ! tidal harmonics data 
    289       INTEGER, OPTIONAL, INTENT(in   ) ::   kit         ! Barotropic timestep counter (for timesplitting option) 
    290       INTEGER, OPTIONAL, INTENT(in   ) ::   kt_offset   ! time offset in units of timesteps. NB. if kit 
    291       !                                                 ! is present then units = subcycle timesteps. 
    292       !                                                 ! kt_offset = 0  => get data at "now"    time level 
    293       !                                                 ! kt_offset = -1 => get data at "before" time level 
    294       !                                                 ! kt_offset = +1 => get data at "after"  time level 
    295       !                                                 ! etc. 
    296       ! 
    297       INTEGER  ::   itide, igrd, ib       ! dummy loop indices 
    298       INTEGER  ::   time_add              ! time offset in units of timesteps 
    299       INTEGER, DIMENSION(3) ::   ilen0    ! length of boundary data (from OBC arrays) 
    300       REAL(wp) ::   z_arg, z_sarg, zflag, zramp   ! local scalars     
    301       REAL(wp), DIMENSION(jpmax_harmo) :: z_sist, z_cost 
    302       !!---------------------------------------------------------------------- 
    303       ! 
    304       ilen0(1) =  SIZE(td%ssh(:,1,1)) 
    305       ilen0(2) =  SIZE(td%u(:,1,1)) 
    306       ilen0(3) =  SIZE(td%v(:,1,1)) 
    307  
    308       zflag=1 
    309       IF ( PRESENT(kit) ) THEN 
    310         IF ( kit /= 1 ) zflag=0 
    311       ENDIF 
    312  
    313       IF ( (nsec_day == NINT(0.5_wp * rdt) .OR. kt==nit000) .AND. zflag==1 ) THEN 
    314         ! 
    315         kt_tide = kt - (nsec_day - 0.5_wp * rdt)/rdt 
    316         ! 
    317         IF(lwp) THEN 
    318            WRITE(numout,*) 
    319            WRITE(numout,*) 'bdytide_update : (re)Initialization of the tidal bdy forcing at kt=',kt 
    320            WRITE(numout,*) '~~~~~~~~~~~~~~ ' 
    321         ENDIF 
    322         ! 
    323         CALL tide_init_elevation ( idx, td ) 
    324         CALL tide_init_velocities( idx, td ) 
    325         ! 
    326       ENDIF  
    327  
    328       time_add = 0 
    329       IF( PRESENT(kt_offset) ) THEN 
    330          time_add = kt_offset 
    331       ENDIF 
    332           
    333       IF( PRESENT(kit) ) THEN   
    334          z_arg = ((kt-kt_tide) * rdt + (kit+0.5_wp*(time_add-1)) * rdt / REAL(nn_baro,wp) ) 
    335       ELSE                               
    336          z_arg = ((kt-kt_tide)+time_add) * rdt 
    337       ENDIF 
    338  
    339       ! Linear ramp on tidal component at open boundaries  
    340       zramp = 1._wp 
    341       IF (ln_tide_ramp) zramp = MIN(MAX( (z_arg + (kt_tide-nit000)*rdt)/(rdttideramp*rday),0._wp),1._wp) 
    342  
    343       DO itide = 1, nb_harmo 
    344          z_sarg = z_arg * omega_tide(itide) 
    345          z_cost(itide) = COS( z_sarg ) 
    346          z_sist(itide) = SIN( z_sarg ) 
    347       END DO 
    348  
    349       DO itide = 1, nb_harmo 
    350          igrd=1                              ! SSH on tracer grid 
    351          DO ib = 1, ilen0(igrd) 
    352             dta%ssh(ib) = dta%ssh(ib) + zramp*(td%ssh(ib,itide,1)*z_cost(itide) + td%ssh(ib,itide,2)*z_sist(itide)) 
    353          END DO 
    354          igrd=2                              ! U grid 
    355          DO ib = 1, ilen0(igrd) 
    356             dta%u2d(ib) = dta%u2d(ib) + zramp*(td%u  (ib,itide,1)*z_cost(itide) + td%u  (ib,itide,2)*z_sist(itide)) 
    357          END DO 
    358          igrd=3                              ! V grid 
    359          DO ib = 1, ilen0(igrd)  
    360             dta%v2d(ib) = dta%v2d(ib) + zramp*(td%v  (ib,itide,1)*z_cost(itide) + td%v  (ib,itide,2)*z_sist(itide)) 
    361          END DO 
    362       END DO 
    363       ! 
    364    END SUBROUTINE bdytide_update 
    365269 
    366270 
     
    406310      ! Linear ramp on tidal component at open boundaries  
    407311      zramp = 1. 
    408       IF (ln_tide_ramp) zramp = MIN(MAX( (z_arg - nit000*rdt)/(rdttideramp*rday),0.),1.) 
     312      IF (ln_tide_ramp) zramp = MIN(MAX( (z_arg - nit000*rdt)/(rn_tide_ramp_dt*rday),0.),1.) 
    409313 
    410314      DO ib_bdy = 1,nb_bdy 
     
    447351            DO itide = 1, nb_harmo 
    448352               ! 
    449                z_sarg = (z_arg + zoff) * omega_tide(itide) 
     353               z_sarg = (z_arg + zoff) * tide_harmonics(itide)%omega 
    450354               z_cost = zramp * COS( z_sarg ) 
    451355               z_sist = zramp * SIN( z_sarg ) 
     
    505409         END DO 
    506410         DO ib = 1 , ilen0(igrd) 
    507             mod_tide(ib)=mod_tide(ib)*ftide(itide) 
    508             phi_tide(ib)=phi_tide(ib)+v0tide(itide)+utide(itide) 
     411            mod_tide(ib)=mod_tide(ib)*tide_harmonics(itide)%f 
     412            phi_tide(ib)=phi_tide(ib)+tide_harmonics(itide)%v0+tide_harmonics(itide)%u 
    509413         ENDDO 
    510414         DO ib = 1 , ilen0(igrd) 
     
    544448         END DO 
    545449         DO ib = 1, ilen0(igrd) 
    546             mod_tide(ib)=mod_tide(ib)*ftide(itide) 
    547             phi_tide(ib)=phi_tide(ib)+v0tide(itide)+utide(itide) 
     450            mod_tide(ib)=mod_tide(ib)*tide_harmonics(itide)%f 
     451            phi_tide(ib)=phi_tide(ib)+tide_harmonics(itide)%v0 + tide_harmonics(itide)%u 
    548452         ENDDO 
    549453         DO ib = 1, ilen0(igrd) 
     
    565469         END DO 
    566470         DO ib = 1, ilen0(igrd) 
    567             mod_tide(ib)=mod_tide(ib)*ftide(itide) 
    568             phi_tide(ib)=phi_tide(ib)+v0tide(itide)+utide(itide) 
     471            mod_tide(ib)=mod_tide(ib)*tide_harmonics(itide)%f 
     472            phi_tide(ib)=phi_tide(ib)+tide_harmonics(itide)%v0 + tide_harmonics(itide)%u 
    569473         ENDDO 
    570474         DO ib = 1, ilen0(igrd) 
  • NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/DYN/dynspg.F90

    r11960 r12205  
    2121   USE dynspg_exp     ! surface pressure gradient     (dyn_spg_exp routine) 
    2222   USE dynspg_ts      ! surface pressure gradient     (dyn_spg_ts  routine) 
    23    USE sbctide        !  
    24    USE updtide        !  
     23   USE tide_mod       ! 
    2524   USE trd_oce        ! trends: ocean variables 
    2625   USE trddyn         ! trend manager: dynamics 
     
    4342   INTEGER, PARAMETER ::   np_EXP = 0   !       explicit time stepping 
    4443   INTEGER, PARAMETER ::   np_NO  =-1   ! no surface pressure gradient, no scheme 
     44   ! 
     45   REAL(wp) ::   zt0step !   Time of day at the beginning of the time step 
    4546 
    4647   !! * Substitutions 
     
    116117         IF( .NOT.ln_dynspg_ts .AND. ( ln_tide_pot .AND. ln_tide )  ) THEN   ! N.B. added directly at sub-time-step in ts-case 
    117118            ! 
    118             CALL upd_tide( kt )                      ! update tide potential 
     119            ! Update tide potential at the beginning of current time step 
     120            zt0step = REAL(nsec_day, wp)-0.5_wp*rdt 
     121            CALL upd_tide(zt0step, Kmm) 
    119122            ! 
    120123            DO jj = 2, jpjm1                         ! add tide potential forcing 
  • NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/DYN/dynspg_ts.F90

    r12150 r12205  
    4444   USE bdytides        ! open boundary condition data 
    4545   USE bdydyn2d        ! open boundary conditions on barotropic variables 
    46    USE sbctide         ! tides 
    47    USE updtide         ! tide potential 
     46   USE tide_mod        ! 
    4847   USE sbcwave         ! surface wave 
    4948   USE diatmb          ! Top,middle,bottom output 
     
    173172      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: ztwdmask, zuwdmask, zvwdmask ! ROMS wetting and drying masks at t,u,v points 
    174173      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zuwdav2, zvwdav2    ! averages over the sub-steps of zuwdmask and zvwdmask 
     174      REAL(wp) ::   zt0substep !   Time of day at the beginning of the time substep 
    175175      !!---------------------------------------------------------------------- 
    176176      ! 
     
    444444         ! 
    445445         IF( ln_bdy      .AND. ln_tide )   CALL bdy_dta_tides( kt, kit=jn, kt_offset= noffset+1 ) 
    446          IF( ln_tide_pot .AND. ln_tide )   CALL upd_tide     ( kt, kit=jn, kt_offset= noffset   ) 
     446         ! Update tide potential at the beginning of current time substep 
     447         IF( ln_tide_pot .AND. ln_tide ) THEN 
     448            zt0substep = REAL(nsec_day, wp) - 0.5_wp*rdt + (jn + noffset - 1) * rdt / REAL(nn_baro, wp) 
     449            CALL upd_tide(zt0substep, Kmm) 
     450         END IF 
    447451         ! 
    448452         !                    !==  extrapolation at mid-step  ==!   (jn+1/2) 
  • NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/IOM/iom.F90

    r12193 r12205  
    5555   LOGICAL, PUBLIC, PARAMETER ::   lk_iomput = .FALSE.       !: iom_put flag 
    5656#endif 
    57    PUBLIC iom_init, iom_swap, iom_open, iom_close, iom_setkt, iom_varid, iom_get, iom_get_var 
     57   PUBLIC iom_init, iom_init_closedef, iom_swap, iom_open, iom_close, iom_setkt, iom_varid, iom_get, iom_get_var 
    5858   PUBLIC iom_chkatt, iom_getatt, iom_putatt, iom_getszuld, iom_rstput, iom_delay_rst, iom_put 
    59    PUBLIC iom_use, iom_context_finalize, iom_miss_val 
     59   PUBLIC iom_use, iom_context_finalize, iom_update_file_name, iom_miss_val 
    6060 
    6161   PRIVATE iom_rp0d, iom_rp1d, iom_rp2d, iom_rp3d 
     
    6464#if defined key_iomput 
    6565   PRIVATE iom_set_domain_attr, iom_set_axis_attr, iom_set_field_attr, iom_set_file_attr, iom_get_file_attr, iom_set_grid_attr 
    66    PRIVATE set_grid, set_grid_bounds, set_scalar, set_xmlatt, set_mooring, iom_update_file_name, iom_sdate 
     66   PRIVATE set_grid, set_grid_bounds, set_scalar, set_xmlatt, set_mooring, iom_sdate 
    6767   PRIVATE iom_set_rst_context, iom_set_rstw_active, iom_set_rstr_active 
    6868# endif 
     
    9292CONTAINS 
    9393 
    94    SUBROUTINE iom_init( cdname, fname, ld_tmppatch )  
     94   SUBROUTINE iom_init( cdname, fname, ld_tmppatch, ld_closedef )  
    9595      !!---------------------------------------------------------------------- 
    9696      !!                     ***  ROUTINE   *** 
     
    102102      CHARACTER(len=*), OPTIONAL, INTENT(in)  :: fname 
    103103      LOGICAL         , OPTIONAL, INTENT(in)  :: ld_tmppatch 
     104      LOGICAL         , OPTIONAL, INTENT(in)  :: ld_closedef 
    104105#if defined key_iomput 
    105106      ! 
     
    116117      INTEGER ::   nldi_save, nlei_save    !:      and close boundaries in output files 
    117118      INTEGER ::   nldj_save, nlej_save    !: 
     119      LOGICAL ::   ll_closedef = .TRUE. 
    118120      !!---------------------------------------------------------------------- 
    119121      ! 
     
    130132         IF( njmpp + jpj - 1 == jpjglo ) nlej = jpj 
    131133      ENDIF 
     134      IF ( PRESENT(ld_closedef) ) ll_closedef = ld_closedef 
    132135      ! 
    133136      ALLOCATE( zt_bnds(2,jpk), zw_bnds(2,jpk) ) 
     
    268271      ENDIF 
    269272      ! 
    270       ! end file definition 
     273      ! set time step length 
    271274      dtime%second = rdt 
    272275      CALL xios_set_timestep( dtime ) 
    273       CALL xios_close_context_definition() 
    274       CALL xios_update_calendar( 0 ) 
     276      ! 
     277      ! conditional closure of context definition 
     278      IF ( ll_closedef ) CALL iom_init_closedef 
    275279      ! 
    276280      DEALLOCATE( zt_bnds, zw_bnds ) 
     
    283287      ! 
    284288   END SUBROUTINE iom_init 
     289 
     290   SUBROUTINE iom_init_closedef 
     291      !!---------------------------------------------------------------------- 
     292      !!            ***  SUBROUTINE iom_init_closedef  *** 
     293      !!---------------------------------------------------------------------- 
     294      !! 
     295      !! ** Purpose : Closure of context definition 
     296      !! 
     297      !!---------------------------------------------------------------------- 
     298 
     299      CALL xios_close_context_definition() 
     300      CALL xios_update_calendar( 0 ) 
     301 
     302   END SUBROUTINE iom_init_closedef 
    285303 
    286304   SUBROUTINE iom_set_rstw_var_active(field) 
  • NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/TDE/tide_mod.F90

    r12065 r12205  
    1717   !!---------------------------------------------------------------------- 
    1818 
    19    USE oce, ONLY : sshn       ! sea-surface height 
     19   USE oce, ONLY : ssh        ! sea-surface height 
    2020   USE par_oce                ! ocean parameters 
    2121   USE phycst, ONLY : rpi     ! pi 
     
    132132      sn_tide_cnames(:) = '' 
    133133      ! Read Namelist nam_tide 
    134       REWIND( numnam_ref )              ! Namelist nam_tide in reference namelist : Tides 
    135134      READ  ( numnam_ref, nam_tide, IOSTAT = ios, ERR = 901) 
    136135901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'nam_tide in reference namelist' ) 
    137136      ! 
    138       REWIND( numnam_cfg )              ! Namelist nam_tide in configuration namelist : Tides 
    139137      READ  ( numnam_cfg, nam_tide, IOSTAT = ios, ERR = 902 ) 
    140138902   IF( ios >  0 )   CALL ctl_nam ( ios , 'nam_tide in configuration namelist' ) 
     
    715713 
    716714 
    717    SUBROUTINE upd_tide(pdelta) 
     715   SUBROUTINE upd_tide(pdelta, Kmm) 
    718716      !!---------------------------------------------------------------------- 
    719717      !!                 ***  ROUTINE upd_tide  *** 
     
    726724      !!----------------------------------------------------------------------       
    727725      REAL, INTENT(in)              ::   pdelta      ! Temporal offset in seconds 
     726      INTEGER, INTENT(IN)           ::   Kmm         ! Time level index 
    728727      INTEGER                       ::   jk          ! Dummy loop index 
    729728      REAL(wp)                      ::   zt, zramp   ! Local scalars 
     
    755754      ! 
    756755      IF( ln_tide_dia ) THEN          ! Output total tidal potential (incl. load potential) 
    757          IF ( iom_use( "tide_pot" ) ) CALL iom_put( "tide_pot", pot_astro(:,:) + rn_scal_load * sshn(:,:) ) 
     756         IF ( iom_use( "tide_pot" ) ) CALL iom_put( "tide_pot", pot_astro(:,:) + rn_scal_load * ssh(:,:,Kmm) ) 
    758757      END IF 
    759758      ! 
  • NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/nemogcm.F90

    r12193 r12205  
    4646   USE closea         ! treatment of closed seas (for ln_closea) 
    4747   USE usrdef_nam     ! user defined configuration 
    48    USE tideini        ! tidal components initialization   (tide_ini routine) 
     48   USE tide_mod, ONLY : tide_init ! tidal components initialization   (tide_init routine) 
    4949   USE bdy_oce,  ONLY : ln_bdy 
    5050   USE bdyini         ! open boundary cond. setting       (bdy_init routine) 
     
    5959   USE diaobs         ! Observation diagnostics       (dia_obs_init routine) 
    6060   USE diacfl         ! CFL diagnostics               (dia_cfl_init routine) 
    61    USE diaharm        ! tidal harmonics diagnostics  (dia_harm_init routine) 
     61   USE diamlr         ! IOM context management for multiple-linear-regression analysis 
    6262   USE step           ! NEMO time-stepping                 (stp     routine) 
    6363   USE isfstp         ! ice shelf                     (isf_stp_init routine) 
     
    7575   USE diatmb         ! Top,middle,bottom output 
    7676   USE dia25h         ! 25h mean output 
     77   USE diadetide      ! Weights computation for daily detiding of model diagnostics 
    7778   USE sbc_oce , ONLY : lk_oasis 
    7879   USE wet_dry        ! Wetting and drying setting   (wad_init routine) 
     
    492493                           CALL dia_tmb_init    ! TMB outputs 
    493494                           CALL dia_25h_init( Nbb )    ! 25h mean  outputs 
    494                            CALL dia_harm_init   ! tidal harmonics outputs 
     495                           CALL dia_detide_init ! Weights computation for daily detiding of model diagnostics 
    495496      IF( ln_diaobs    )   CALL dia_obs( nit000-1, Nnn )   ! Observation operator for restart 
     497                           CALL dia_mlr_init    ! Initialisation of IOM context management for multiple-linear-regression analysis 
    496498 
    497499      !                                      ! Assimilation increments 
  • NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/step.F90

    r12193 r12205  
    107107                              
    108108      IF( kstp == nit000 ) THEN                       ! initialize IOM context (must be done after nemo_init for AGRIF+XIOS+OASIS) 
    109                              CALL iom_init(      cxios_context          )  ! for model grid (including passible AGRIF zoom) 
     109                             CALL iom_init( cxios_context, ld_closedef=.FALSE. )   ! for model grid (including passible AGRIF zoom) 
     110         IF( lk_diamlr   )   CALL dia_mlr_iom_init    ! with additional setup for multiple-linear-regression analysis 
     111                             CALL iom_init_closedef 
    110112         IF( ln_crs      )   CALL iom_init( TRIM(cxios_context)//"_crs" )  ! for coarse grid 
    111113      ENDIF 
     
    117119      ! Update external forcing (tides, open boundaries, ice shelf interaction and surface boundary condition (including sea-ice) 
    118120      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    119       IF( ln_tide    )   CALL sbc_tide( kstp )                        ! update tide potential 
     121      IF( ln_tide    )   CALL tide_update( kstp )                     ! update tide potential 
    120122      IF( ln_apr_dyn )   CALL sbc_apr ( kstp )                        ! atmospheric pressure (NB: call before bdy_dta which needs ssh_ib)  
    121123      IF( ln_bdy     )   CALL bdy_dta ( kstp, Nnn, kt_offset = +1 )   ! update dynamic & tracer data at open boundaries 
     
    211213      ! diagnostics and outputs 
    212214      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    213       IF( ln_floats  )   CALL flo_stp ( kstp, Nbb, Nnn )        ! drifting Floats 
    214       IF( ln_diacfl  )   CALL dia_cfl ( kstp,      Nnn )        ! Courant number diagnostics 
    215                          CALL dia_hth ( kstp,      Nnn )        ! Thermocline depth (20 degres isotherm depth) 
    216       IF( ln_diadct  )   CALL dia_dct ( kstp,      Nnn )        ! Transports 
    217                          CALL dia_ar5 ( kstp,      Nnn )        ! ar5 diag 
    218                          CALL dia_ptr ( kstp,      Nnn )        ! Poleward adv/ldf TRansports diagnostics 
    219       IF( ln_diaharm )   CALL dia_harm( kstp,      Nnn )        ! Tidal harmonic analysis 
    220                          CALL dia_wri ( kstp,      Nnn )        ! ocean model: outputs 
    221       ! 
    222       IF( ln_crs     )   CALL crs_fld ( kstp,      Nnn )   ! ocean model: online field coarsening & output 
     215      IF( ln_floats  )   CALL flo_stp   ( kstp, Nbb, Nnn )      ! drifting Floats 
     216      IF( ln_diacfl  )   CALL dia_cfl   ( kstp,      Nnn )      ! Courant number diagnostics 
     217                         CALL dia_hth   ( kstp,      Nnn )      ! Thermocline depth (20 degres isotherm depth) 
     218      IF( ln_diadct  )   CALL dia_dct   ( kstp,      Nnn )      ! Transports 
     219                         CALL dia_ar5   ( kstp,      Nnn )      ! ar5 diag 
     220                         CALL dia_ptr   ( kstp,      Nnn )      ! Poleward adv/ldf TRansports diagnostics 
     221                         CALL dia_wri   ( kstp,      Nnn )      ! ocean model: outputs 
     222      IF( ln_crs     )   CALL crs_fld   ( kstp,      Nnn )      ! ocean model: online field coarsening & output 
     223      IF( lk_diadetide ) CALL dia_detide( kstp )                ! Weights computation for daily detiding of model diagnostics 
     224      IF( lk_diamlr  )   CALL dia_mlr                           ! Update time used in multiple-linear-regression analysis 
    223225       
    224226#if defined key_top 
  • NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/step_oce.F90

    r12150 r12205  
    1919   USE sbccpl          ! surface boundary condition: coupled formulation (call send at end of step) 
    2020   USE sbcapr          ! surface boundary condition: atmospheric pressure 
    21    USE sbctide         ! Tide initialisation 
     21   USE tide_mod, ONLY : ln_tide, tide_update 
    2222   USE sbcwave         ! Wave intialisation 
    2323 
     
    8181   USE diahth          ! thermocline depth                (dia_hth routine) 
    8282   USE diahsb          ! heat, salt and volume budgets    (dia_hsb routine) 
    83    USE diaharm 
    8483   USE diacfl 
    8584   USE diaobs          ! Observation operator 
     85   USE diadetide       ! Weights computation for daily detiding of model diagnostics 
     86   USE diamlr          ! IOM context management for multiple-linear-regression analysis 
    8687   USE flo_oce         ! floats variables 
    8788   USE floats          ! floats computation               (flo_stp routine) 
Note: See TracChangeset for help on using the changeset viewer.