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 9403 for branches/2017/dev_merge_2017/NEMOGCM/CONFIG/TEST_CASES/CANAL/MY_SRC/usrdef_istate.F90 – NEMO

Ignore:
Timestamp:
2018-03-15T15:57:42+01:00 (6 years ago)
Author:
smasson
Message:

dev_merge_2017: update test_cases CANAL

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2017/dev_merge_2017/NEMOGCM/CONFIG/TEST_CASES/CANAL/MY_SRC/usrdef_istate.F90

    r9302 r9403  
    33   !!                     ***  MODULE usrdef_istate   *** 
    44   !! 
    5    !!                      ===  VORTEX configuration  === 
     5   !!                      ===  CANAL configuration  === 
    66   !! 
    77   !! User defined : set the initial state of a user configuration 
     
    1414   !!---------------------------------------------------------------------- 
    1515   USE par_oce        ! ocean space and time domain 
    16    USE dom_oce , ONLY : glamt, gphit, glamu, gphiu, glamv, gphiv, ff_t, ff_f 
     16   USE dom_oce         
    1717   USE phycst         ! physical constants 
    1818   ! 
     
    2121   USE lbclnk         ! ocean lateral boundary conditions (or mpp link) 
    2222   !    
    23    USE usrdef_nam, ONLY : rn_ppgphi0, rn_vtxmax, rn_uzonal, rn_ujetszx, rn_ujetszy, nn_initcase, ln_sshnoise !  
     23   USE usrdef_nam 
    2424    
    2525   IMPLICIT NONE 
     
    4040      !!  
    4141      !! ** Purpose :   Initialization of the dynamics and tracers 
    42       !!                Here VORTEX configuration  
     42      !!                Here CANAL configuration  
    4343      !! 
    4444      !! ** Method  :   Set a gaussian anomaly of pressure and associated 
     
    5656      REAL(wp) :: zpsurf, zdyPs, zdxPs 
    5757      REAL(wp) :: zdt, zdu, zdv 
    58       REAL(wp) :: zjetx, zjety 
     58      REAL(wp) :: zjetx, zjety, zbeta 
    5959      REAL(wp), DIMENSION(jpi,jpj)  ::   zrandom 
    6060      !!---------------------------------------------------------------------- 
    6161      ! 
    6262      IF(lwp) WRITE(numout,*) 
    63       IF(lwp) WRITE(numout,*) 'usr_def_istate : VORTEX configuration, analytical definition of initial state' 
     63      IF(lwp) WRITE(numout,*) 'usr_def_istate : CANAL configuration, analytical definition of initial state' 
    6464      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~   ' 
    6565      ! 
     
    7474         pssh(:,:) = 0. 
    7575         ! temperature: 
    76          pts(:,:,:,jp_tem) = 10._wp * ptmask(:,:,:)  
    77          ! salinity:   
    78          pts(:,:,:,jp_sal) = 35._wp * ptmask(:,:,:)  
     76         pts(:,:,:,jp_tem) = 10._wp 
     77         ! salinity:   
     78         pts(:,:,:,jp_sal) = 35._wp 
    7979         ! velocities: 
    8080         pu(:,:,:) = 0. 
     
    8282          
    8383      CASE(1)    ! geostrophic zonal jet from -zjety to +zjety 
    84           
    85          ! sea level: ssh = - fuy / g 
    86          WHERE( ABS(gphit) <= zjety ) 
    87             pssh(:,:) = - ff_t(:,:) * rn_uzonal * gphit(:,:) * 1.e3 / grav  * ptmask(:,:,1) 
    88          ELSEWHERE 
    89             pssh(:,:) = - ff_t(:,:) * rn_uzonal * SIGN(zjety, gphit(:,:)) * 1.e3 / grav  * ptmask(:,:,1) 
    90          END WHERE 
    91          ! temperature: 
    92          pts(:,:,:,jp_tem) = 10._wp * ptmask(:,:,:)  
    93          ! salinity:   
    94          pts(:,:,:,jp_sal) = 1._wp * ptmask(:,:,:)  
    95          DO jk=1, jpkm1 
    96             WHERE( ABS(gphit) <= zjety ) pts(:,:,jk,jp_sal) = 2. 
     84 
     85         ! sea level: 
     86         SELECT CASE( nn_fcase ) 
     87         CASE(0)    ! f = f0 
     88            ! sea level: ssh = - fuy / g 
     89            WHERE( ABS(gphit) <= zjety ) 
     90               pssh(:,:) = - ff_t(:,:) * rn_uzonal * gphit(:,:) * 1.e3 / grav 
     91            ELSEWHERE 
     92               pssh(:,:) = - ff_t(:,:) * rn_uzonal * SIGN(zjety, gphit(:,:)) * 1.e3 / grav 
     93            END WHERE 
     94         CASE(1)    ! f = f0 + beta*y 
     95            ! sea level: ssh = - u / g * ( fy + 0.5 * beta * y^2 ) 
     96            zbeta = 2._wp * omega * COS( rad * rn_ppgphi0 ) / ra 
     97            WHERE( ABS(gphit) <= zjety ) 
     98               pssh(:,:) = - rn_uzonal / grav * ( ff_t(:,:) * gphit(:,:) * 1.e3 + 0.5 * zbeta * gphit(:,:) * gphit(:,:) * 1.e6 ) 
     99            ELSEWHERE 
     100               pssh(:,:) = - rn_uzonal / grav * ( ff_t(:,:) * SIGN(zjety, gphit(:,:)) * 1.e3   & 
     101                  &                             + 0.5 * zbeta * zjety * zjety * 1.e6 ) 
     102            END WHERE 
     103         END SELECT 
     104         ! temperature: 
     105         pts(:,:,:,jp_tem) = 10._wp 
     106         ! salinity:   
     107         pts(:,:,jpk,jp_sal) = 0. 
     108         DO jk=1, jpkm1 
     109            pts(:,:,jk,jp_sal) = gphit(:,:) 
    97110         END DO 
    98111         ! velocities: 
     
    102115         END DO 
    103116         pv(:,:,:) = 0. 
    104           
    105       CASE(10)    ! geostrophic zonal current shear 
    106           
    107          ! sea level: ssh = - fuy / g 
    108          WHERE( ABS(gphit) <= zjety ) 
    109             pssh(:,:) = - ff_t(:,:) * rn_uzonal * ABS(gphit(:,:) * 1.e3) / grav  * ptmask(:,:,1) 
    110          ELSEWHERE 
    111             pssh(:,:) = - ff_t(:,:) * rn_uzonal * zjety * 1.e3 / grav  * ptmask(:,:,1) 
    112          END WHERE 
    113          ! temperature: 
    114          pts(:,:,:,jp_tem) = 10._wp * ptmask(:,:,:)  
    115          ! salinity:   
    116          pts(:,:,:,jp_sal) = 2._wp * ptmask(:,:,:)  
    117          DO jk=1, jpkm1 
    118             WHERE( ABS(gphit) <= zjety ) pts(:,:,jk,jp_sal) = 2. + SIGN(1.,gphit(:,:)) 
     117         !                   
     118      CASE(2)    ! geostrophic zonal current shear 
     119       
     120         ! sea level: 
     121         SELECT CASE( nn_fcase ) 
     122         CASE(0)    ! f = f0 
     123            ! sea level: ssh = - fuy / g 
     124            WHERE( ABS(gphit) <= zjety ) 
     125               pssh(:,:) = - ff_t(:,:) * rn_uzonal * ABS(gphit(:,:)) * 1.e3 / grav 
     126            ELSEWHERE 
     127               pssh(:,:) = - ff_t(:,:) * rn_uzonal * zjety * 1.e3 / grav 
     128            END WHERE 
     129         CASE(1)    ! f = f0 + beta*y 
     130            ! sea level: ssh = - u / g * ( fy + 0.5 * beta * y^2 ) 
     131            zbeta = 2._wp * omega * COS( rad * rn_ppgphi0 ) / ra 
     132            WHERE( ABS(gphit) <= zjety ) 
     133               pssh(:,:) = - SIGN(rn_uzonal, gphit(:,:)) / grav   & 
     134                  &        * ( ff_t(:,:) * gphit(:,:) * 1.e3 + 0.5 * zbeta * gphit(:,:) * gphit(:,:) * 1.e6 ) 
     135            ELSEWHERE 
     136               pssh(:,:) = - SIGN(rn_uzonal, gphit(:,:)) / grav   & 
     137                  &        * ( ff_t(:,:) * SIGN(zjety, gphit(:,:)) * 1.e3 + 0.5 * zbeta * zjety * zjety * 1.e6 ) 
     138            END WHERE 
     139         END SELECT 
     140         ! temperature: 
     141         pts(:,:,:,jp_tem) = 10._wp 
     142         ! salinity:   
     143         pts(:,:,:,jp_sal) = 2. 
     144         DO jk=1, jpkm1 
     145            WHERE( ABS(gphiv) <= zjety ) pts(:,:,jk,jp_sal) = 2. + SIGN(1.,gphiv(:,:)) 
    119146         END DO 
    120147         ! velocities: 
    121148         pu(:,:,:) = 0. 
    122149         DO jk=1, jpkm1 
    123             WHERE( ABS(gphit) <= zjety ) pu(:,:,jk) = SIGN(rn_uzonal,gphit(:,:))*SIGN(1.,rn_uzonal) 
    124             WHERE( ABS(gphit) == 0.    ) pu(:,:,jk) = 0.   
    125          END DO 
    126          pv(:,:,:) = 0. 
    127           
    128       CASE(11)    ! geostrophic zonal pulse 
     150            WHERE( ABS(gphiv) <= zjety ) pu(:,:,jk) = SIGN(rn_uzonal,gphit(:,:))*SIGN(1.,rn_uzonal) 
     151            WHERE( ABS(gphiv) == 0.    ) pu(:,:,jk) = 0.   
     152         END DO 
     153         pv(:,:,:) = 0. 
     154         !                   
     155      CASE(3)    ! gaussian zonal currant 
     156 
     157         ! zonal current 
     158         DO jk=1, jpkm1 
     159            ! gphit and lambda are both in km 
     160            pu(:,:,jk) = rn_uzonal * EXP( - 0.5 * gphit(:,:)**2 / rn_lambda**2 ) 
     161         END DO 
     162          
     163         ! sea level: 
     164         pssh(:,1) = - ff_t(:,1) / grav * pu(:,1,1) * e2t(:,1) 
     165         DO jl=1, jpnj 
     166            DO jj=nldj, nlej 
     167               DO ji=nldi, nlei 
     168                  pssh(ji,jj) = pssh(ji,jj-1) - ff_t(ji,jj) / grav * pu(ji,jj,1) * e2t(ji,jj) 
     169               END DO 
     170            END DO 
     171            CALL lbc_lnk( pssh, 'T',  1. ) 
     172         END DO 
     173          
     174         ! temperature: 
     175         pts(:,:,:,jp_tem) = 10._wp 
     176         ! salinity:   
     177         DO jk=1, jpkm1 
     178            pts(:,:,jk,jp_sal) = gphit(:,:) 
     179         END DO 
     180         ! velocities: 
     181         pv(:,:,:) = 0. 
     182         !             
     183      CASE(4)    ! geostrophic zonal pulse 
    129184    
    130185         DO jj=1, jpj 
     
    138193               END IF 
    139194               IF ( ABS(gphit(ji,jj)) <= zjety ) THEN 
    140                   pssh(ji,jj) = - ff_t(ji,jj) * zdu * gphit(ji,jj) * 1.e3 / grav  * ptmask(ji,jj,1) 
     195                  pssh(ji,jj) = - ff_t(ji,jj) * zdu * gphit(ji,jj) * 1.e3 / grav 
    141196                  pu(ji,jj,:) = zdu 
    142197                  pts(ji,jj,:,jp_sal) = zdu / rn_uzonal + 1. 
    143198               ELSE 
    144                   pssh(ji,jj) = - ff_t(ji,jj) * zdu * SIGN(zjety,gphit(ji,jj)) * 1.e3 / grav  * ptmask(ji,jj,1) 
     199                  pssh(ji,jj) = - ff_t(ji,jj) * zdu * SIGN(zjety,gphit(ji,jj)) * 1.e3 / grav  
    145200                  pu(ji,jj,:) = 0. 
    146201                  pts(ji,jj,:,jp_sal) = 1. 
     
    153208         pv(:,:,:) = 0. 
    154209          
    155        CASE(2)    ! vortex 
     210          
     211       CASE(5)    ! vortex 
    156212                  ! 
    157213         zf0   = 2._wp * omega * SIN( rad * rn_ppgphi0 ) 
    158214         zumax = rn_vtxmax * SIGN(1._wp, zf0) ! Here Anticyclonic: set zumax=-1 for cyclonic 
    159          zlambda = SQRT(2._wp)*60.e3        ! Horizontal scale in meters  
     215         zlambda = SQRT(2._wp)*rn_lambda       ! Horizontal scale in meters  
    160216         zn2 = 3.e-3**2 
    161217         zH = 0.5_wp * 5000._wp 
     
    233289      END SELECT 
    234290 
    235       IF (ln_sshnoise) pssh(:,:) = pssh(:,:) + (0.1*zrandom-0.05) 
     291      IF (ln_sshnoise) CALL RANDOM_NUMBER(zrandom) 
     292      IF (ln_sshnoise) THEN 
     293         CALL RANDOM_NUMBER(zrandom) 
     294         pssh(:,:) = pssh(:,:) + ( 0.1  * zrandom(:,:) - 0.05 ) 
     295      END IF 
    236296      CALL lbc_lnk( pssh, 'T',  1. ) 
    237297      CALL lbc_lnk(  pts, 'T',  1. ) 
Note: See TracChangeset for help on using the changeset viewer.