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 11268 – NEMO

Changeset 11268


Ignore:
Timestamp:
2019-07-15T15:20:49+02:00 (5 years ago)
Author:
smasson
Message:

dev_r10984_HPC-13 : introduce a logical to force vertical interpolation if the number of vertical levels in bdy dta == jpk, see #2285

Location:
NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization
Files:
6 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/cfgs/AMM12/EXPREF/namelist_cfg

    r11267 r11268  
    187187&nambdy_dta    !  open boundaries - external data                        
    188188!----------------------------------------------------------------------- 
    189    ln_full_vel = .false. 
    190  
     189   ln_zinterp  = .false.      !  T if a vertical interpolation is required. Variables gdep[tuv] and e3[tuv] must exist in the file 
     190   !                          !  automatically defined to T if the number of vertical levels in bdy dta /= jpk 
     191   ln_full_vel = .false.      !  T if [uv]3d are "full" velocities and not only its baroclinic components 
     192   !                          !  in this case, baroclinic and barotropic velocities will be recomputed -> [uv]2d not needed 
     193   ! 
    191194   cn_dir  =    './bdydta/' 
    192195   !___________!_________________________!___________________!___________!_____________!________!___________!__________________!__________!_______________! 
  • NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/cfgs/SHARED/namelist_ref

    r11267 r11268  
    605605   ln_dyn3d_dmp  =.false.     !  open boundary condition for baroclinic velocities 
    606606   rn_time_dmp   =  1.        !  Damping time scale in days 
    607    rn_time_dmp_out = 1.        !  Outflow damping time scale 
     607   rn_time_dmp_out = 1.       !  Outflow damping time scale 
    608608   nn_rimwidth   = 10         !  width of the relaxation zone 
    609609   ln_vol        = .false.    !  total volume correction (see nn_volctl parameter) 
     
    613613&nambdy_dta    !  open boundaries - external data                       (see nam_bdy) 
    614614!----------------------------------------------------------------------- 
    615    ln_full_vel = .false.      !  ??? 
    616  
     615   ln_zinterp  = .false.      !  T if a vertical interpolation is required. Variables gdep[tuv] and e3[tuv] must exist in the file 
     616   !                          !  automatically defined to T if the number of vertical levels in bdy dta /= jpk 
     617   ln_full_vel = .false.      !  T if [uv]3d are "full" velocities and not only its baroclinic components 
     618   !                          !  in this case, baroclinic and barotropic velocities will be recomputed -> [uv]2d not needed 
     619   ! 
    617620   cn_dir      = 'bdydta/'    !  root directory for the BDY data location 
    618621   !___________!_________________________!___________________!___________!_____________!________!___________!__________________!__________!_______________! 
  • NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/cfgs/SPITZ12/EXPREF/namelist_cfg

    r11267 r11268  
    181181&nambdy_dta    !  open boundaries - external data                       (see nam_bdy) 
    182182!----------------------------------------------------------------------- 
    183    ln_full_vel = .false.      !  ??? 
    184  
     183   ln_zinterp  = .false.      !  T if a vertical interpolation is required. Variables gdep[tuv] and e3[tuv] must exist in the file 
     184   !                          !  automatically defined to T if the number of vertical levels in bdy dta /= jpk 
     185   ln_full_vel = .false.      !  T if [uv]3d are "full" velocities and not only its baroclinic components 
     186   !                          !  in this case, baroclinic and barotropic velocities will be recomputed -> [uv]2d not needed 
     187   ! 
    185188   cn_dir  =  './' 
    186189!              !  file name      ! frequency (hours) ! variable  ! time interp. !  clim   ! 'yearly'/ ! weights  ! rotation ! land/sea mask ! 
  • NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/BDY/bdydta.F90

    r11224 r11268  
    331331      LOGICAL                                ::   ln_full_vel   ! =T => full velocities in 3D boundary data 
    332332      !                                                         ! =F => baroclinic velocities in 3D boundary data 
    333       INTEGER                                ::   ipk,ipl           ! 
     333      LOGICAL                                ::   ln_zinterp    ! =T => requires a vertical interpolation of the bdydta 
     334      INTEGER                                ::   ipk,ipl       ! 
    334335      INTEGER                                ::   idvar         ! variable ID 
    335336      INTEGER                                ::   indims        ! number of dimensions of the variable 
     
    348349      NAMELIST/nambdy_dta/ cn_dir, bn_tem, bn_sal, bn_u3d, bn_v3d, bn_ssh, bn_u2d, bn_v2d  
    349350      NAMELIST/nambdy_dta/ bn_a_i, bn_h_i, bn_h_s 
    350       NAMELIST/nambdy_dta/ ln_full_vel 
     351      NAMELIST/nambdy_dta/ ln_full_vel, ln_zinterp 
    351352      !!--------------------------------------------------------------------------- 
    352353      ! 
     
    361362      ENDIF 
    362363      bf(:,:)%clrootname = 'NOT USED'   ! default definition used as a flag in fld_read to do nothing. 
    363       bf(:,:)%ltotvel = .FALSE.         ! default definition 
     364      bf(:,:)%lzint      = .FALSE.      ! default definition 
     365      bf(:,:)%ltotvel    = .FALSE.      ! default definition 
    364366  
    365367      ! Read namelists 
     
    521523                  CALL fld_fill( bf_alias, bn_alias, cn_dir, 'bdy_dta', cl3//' '//ctmp1, ctmp2 )   ! use namelist info 
    522524                  IF( bf_alias(1)%ln_tint ) ALLOCATE( bf_alias(1)%fdta( iszdim, 1, ipk, 2 ) ) 
    523                   bf_alias(1)%imap => idx_bdy(jbdy)%nbmap(1:iszdim,igrd)      ! associate the mapping used for this bdy 
    524                   bf_alias(1)%igrd = igrd                                     ! used only for vertical integration of 3D arrays 
     525                  bf_alias(1)%imap    => idx_bdy(jbdy)%nbmap(1:iszdim,igrd)   ! associate the mapping used for this bdy 
     526                  bf_alias(1)%igrd    = igrd                                  ! used only for vertical integration of 3D arrays 
    525527                  bf_alias(1)%ltotvel = ln_full_vel                           ! T if u3d is full velocity 
     528                  bf_alias(1)%lzint   = ln_zinterp                            ! T if it requires a vertical interpolation 
    526529               ENDIF 
    527530 
  • NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/SBC/fldread.F90

    r11267 r11268  
    8686      INTEGER, POINTER, DIMENSION(:)  ::   imap         !   Array of integer pointers to 1D arrays 
    8787      LOGICAL                         ::   ltotvel      !   total velocity or not (T/F) 
     88      LOGICAL                         ::   lzint        !   T if it requires a vertical interpolation 
    8889   END TYPE FLD 
    8990 
     
    616617      ! 
    617618      IF( ASSOCIATED(sdjf%imap) ) THEN 
    618          IF( sdjf%ln_tint ) THEN   ;   CALL fld_map( sdjf%num, sdjf%clvar, sdjf%fdta(:,:,:,2),                & 
    619             &                                        sdjf%nrec_a(1), sdjf%imap, sdjf%igrd, sdjf%ibdy, sdjf%ltotvel ) 
    620          ELSE                      ;   CALL fld_map( sdjf%num, sdjf%clvar, sdjf%fnow(:,:,:  ),                & 
    621             &                                        sdjf%nrec_a(1), sdjf%imap, sdjf%igrd, sdjf%ibdy, sdjf%ltotvel ) 
     619         IF( sdjf%ln_tint ) THEN   ;   CALL fld_map( sdjf%num, sdjf%clvar, sdjf%fdta(:,:,:,2), sdjf%nrec_a(1),   & 
     620            &                                        sdjf%imap, sdjf%igrd, sdjf%ibdy, sdjf%ltotvel, sdjf%lzint ) 
     621         ELSE                      ;   CALL fld_map( sdjf%num, sdjf%clvar, sdjf%fnow(:,:,:  ), sdjf%nrec_a(1),   & 
     622            &                                        sdjf%imap, sdjf%igrd, sdjf%ibdy, sdjf%ltotvel, sdjf%lzint ) 
    622623         ENDIF 
    623624      ELSE IF( LEN(TRIM(sdjf%wgtname)) > 0 ) THEN 
     
    676677 
    677678    
    678    SUBROUTINE fld_map( knum, cdvar, pdta, krec, kmap, kgrd, kbdy, ldtotvel ) 
     679   SUBROUTINE fld_map( knum, cdvar, pdta, krec, kmap, kgrd, kbdy, ldtotvel, ldzint ) 
    679680      !!--------------------------------------------------------------------- 
    680681      !!                    ***  ROUTINE fld_map  *** 
     
    691692      INTEGER, OPTIONAL         , INTENT(in   ) ::   kgrd         ! grid type (t, u, v) 
    692693      INTEGER, OPTIONAL         , INTENT(in   ) ::   kbdy         ! bdy number 
    693       LOGICAL, OPTIONAL         , INTENT(in   ) ::   ldtotvel     ! true if toal ( = barotrop + barocline) velocity 
     694      LOGICAL, OPTIONAL         , INTENT(in   ) ::   ldtotvel     ! true if total ( = barotrop + barocline) velocity 
     695      LOGICAL, OPTIONAL         , INTENT(in   ) ::   ldzint       ! true if 3D variable requires a vertical interpolation 
    694696      !! 
    695697      INTEGER                                   ::   ipi          ! length of boundary data on local process 
     
    707709      CHARACTER(LEN=1),DIMENSION(3)             ::   clgrid 
    708710      LOGICAL                                   ::   lluld        ! is the variable using the unlimited dimension 
     711      LOGICAL                                   ::   llzint       ! local value of ldzint 
    709712      !!--------------------------------------------------------------------- 
    710713      ! 
     
    714717      ipj = SIZE( pdta, 2 )   ! must be equal to 1 
    715718      ipk = SIZE( pdta, 3 ) 
     719      ! 
     720      llzint = .FALSE. 
     721      IF( PRESENT(ldzint) )   llzint = ldzint 
    716722      ! 
    717723      idvar = iom_varid( knum, cdvar, kndims = indims, kdimsz = idimsz, lduld = lluld  ) 
     
    739745         CALL iom_get ( knum, jpdom_unknown, cdvar, zz_read(:,:,:), krec )   ! call iom_get with a 3D file 
    740746         ! 
    741          IF( ipkb /= ipk ) THEN   ! boundary data not on model vertical grid : vertical interpolation 
     747         IF( ipkb /= ipk .OR. llzint ) THEN   ! boundary data not on model vertical grid : vertical interpolation 
    742748            ! 
    743749            IF( ipk == jpk .AND. iom_varid(knum,'gdep'//clgrid(kgrd)) /= -1 .AND. iom_varid(knum,'e3'//clgrid(kgrd)) /= -1 ) THEN 
     
    11791185         sdf(jf)%imap       => NULL() 
    11801186         sdf(jf)%ltotvel    = .FALSE. 
     1187         sdf(jf)%lzint      = .FALSE. 
    11811188      END DO 
    11821189      ! 
  • NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/tests/WAD/EXPREF/namelist_cfg

    r11267 r11268  
    180180&nambdy_dta    !  open boundaries - external data                        
    181181!----------------------------------------------------------------------- 
     182   ln_zinterp  = .false.      !  T if a vertical interpolation is required. Variables gdep[tuv] and e3[tuv] must exist in the file 
     183   !                          !  automatically defined to T if the number of vertical levels in bdy dta /= jpk 
     184   ln_full_vel = .false.      !  T if [uv]3d are "full" velocities and not only its baroclinic components 
     185   !                          !  in this case, baroclinic and barotropic velocities will be recomputed -> [uv]2d not needed 
     186   ! 
     187   cn_dir  =  './' 
    182188!              !  file name      ! frequency (hours) ! variable  ! time interp. !  clim   ! 'yearly'/ ! weights  ! rotation ! land/sea mask ! 
    183189!              !                 !  (if <0  months)  !   name    !  (logical)   !  (T/F ) ! 'monthly' ! filename ! pairing  ! filename      ! 
     
    185191   bn_u2d =     'bdyuv_tc7'  ,         1.       , 'ubdy'  ,     .true.     , .true.  ,  'daily'  ,    ''    ,   ''     , '' 
    186192   bn_v2d =     'bdyuv_tc7'  ,         1.       , 'vbdy'  ,     .true.     , .true.  ,  'daily'  ,    ''    ,   ''     , '' 
    187    cn_dir      =    './'   !  root directory for the location of the bulk files 
    188    ln_full_vel = .false.        ! 
    189193/ 
    190194!----------------------------------------------------------------------- 
Note: See TracChangeset for help on using the changeset viewer.