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.
istate.F90 in branches/2015/dev_r5151_UKMO_ISF/NEMOGCM/NEMO/OPA_SRC/DOM – NEMO

source: branches/2015/dev_r5151_UKMO_ISF/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 5189

Last change on this file since 5189 was 5189, checked in by mathiot, 9 years ago

ISF cleaning branch: simplification and bug correction in hpg_isf, zps_hde_isf, mixed layer definition, slope, diffusion, vertical advection and top friction

  • Property svn:keywords set to Id
File size: 27.4 KB
RevLine 
[3]1MODULE istate
2   !!======================================================================
3   !!                     ***  MODULE  istate  ***
4   !! Ocean state   :  initial state setting
5   !!=====================================================================
[2104]6   !! History :  OPA  !  1989-12  (P. Andrich)  Original code
7   !!            5.0  !  1991-11  (G. Madec)  rewritting
8   !!            6.0  !  1996-01  (G. Madec)  terrain following coordinates
9   !!            8.0  !  2001-09  (M. Levy, M. Ben Jelloul)  istate_eel
10   !!            8.0  !  2001-09  (M. Levy, M. Ben Jelloul)  istate_uvg
11   !!   NEMO     1.0  !  2003-08  (G. Madec, C. Talandier)  F90: Free form, modules + EEL R5
12   !!             -   !  2004-05  (A. Koch-Larrouy)  istate_gyre
13   !!            2.0  !  2006-07  (S. Masson)  distributed restart using iom
14   !!            3.3  !  2010-10  (C. Ethe) merge TRC-TRA
[3294]15   !!            3.4  !  2011-04  (G. Madec) Merge of dtatem and dtasal & suppression of tb,tn/sb,sn
[508]16   !!----------------------------------------------------------------------
[3]17
18   !!----------------------------------------------------------------------
19   !!   istate_init   : initial state setting
20   !!   istate_tem    : analytical profile for initial Temperature
21   !!   istate_sal    : analytical profile for initial Salinity
22   !!   istate_eel    : initial state setting of EEL R5 configuration
[93]23   !!   istate_gyre   : initial state setting of GYRE configuration
[3]24   !!   istate_uvg    : initial velocity in geostropic balance
25   !!----------------------------------------------------------------------
26   USE oce             ! ocean dynamics and active tracers
27   USE dom_oce         ! ocean space and time domain
[4245]28   USE c1d             ! 1D vertical configuration
[2104]29   USE daymod          ! calendar
30   USE eosbn2          ! eq. of state, Brunt Vaisala frequency (eos     routine)
[3]31   USE ldftra_oce      ! ocean active tracers: lateral physics
32   USE zdf_oce         ! ocean vertical physics
33   USE phycst          ! physical constants
[3294]34   USE dtatsd          ! data temperature and salinity   (dta_tsd routine)
[4245]35   USE dtauvd          ! data: U & V current             (dta_uvd routine)
[544]36   USE zpshde          ! partial step: hor. derivative (zps_hde routine)
37   USE eosbn2          ! equation of state            (eos bn2 routine)
[593]38   USE domvvl          ! varying vertical mesh
39   USE dynspg_oce      ! pressure gradient schemes
[4370]40   USE dynspg_flt      ! filtered free surface
[3764]41   USE sol_oce         ! ocean solver variables
[4990]42   !
43   USE in_out_manager  ! I/O manager
44   USE iom             ! I/O library
[2715]45   USE lib_mpp         ! MPP library
[3680]46   USE restart         ! restart
[3294]47   USE wrk_nemo        ! Memory allocation
48   USE timing          ! Timing
[2715]49
[3]50   IMPLICIT NONE
51   PRIVATE
52
[508]53   PUBLIC   istate_init   ! routine called by step.F90
[3]54
55   !! * Substitutions
56#  include "domzgr_substitute.h90"
57#  include "vectopt_loop_substitute.h90"
58   !!----------------------------------------------------------------------
[4990]59   !! NEMO/OPA 3.7 , NEMO Consortium (2014)
[888]60   !! $Id$
[2715]61   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[3]62   !!----------------------------------------------------------------------
63CONTAINS
64
65   SUBROUTINE istate_init
66      !!----------------------------------------------------------------------
67      !!                   ***  ROUTINE istate_init  ***
68      !!
[508]69      !! ** Purpose :   Initialization of the dynamics and tracer fields.
[3]70      !!----------------------------------------------------------------------
[2148]71      ! - ML - needed for initialization of e3t_b
[4354]72      INTEGER  ::  ji,jj,jk     ! dummy loop indices
[4245]73      REAL(wp), POINTER, DIMENSION(:,:,:,:)  ::  zuvd    ! U & V data workspace
[3294]74      !!----------------------------------------------------------------------
75      !
[4990]76      IF( nn_timing == 1 )   CALL timing_start('istate_init')
[3294]77      !
[3]78
[4990]79      IF(lwp) WRITE(numout,*) ' '
[508]80      IF(lwp) WRITE(numout,*) 'istate_ini : Initialization of the dynamics and tracers'
81      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
[3]82
[3294]83      CALL dta_tsd_init                       ! Initialisation of T & S input data
[4245]84      IF( lk_c1d ) CALL dta_uvd_init          ! Initialization of U & V input data
[3]85
[4990]86      rhd  (:,:,:  ) = 0._wp
87      rhop (:,:,:  ) = 0._wp
88      rn2  (:,:,:  ) = 0._wp
89      tsa  (:,:,:,:) = 0._wp   
90      rab_b(:,:,:,:) = 0._wp
91      rab_n(:,:,:,:) = 0._wp
[3294]92
[15]93      IF( ln_rstart ) THEN                    ! Restart from a file
[3]94         !                                    ! -------------------
95         CALL rst_read                           ! Read the restart file
[1130]96         CALL day_init                           ! model calendar (using both namelist and restart infos)
[3]97      ELSE
98         !                                    ! Start from rest
99         !                                    ! ---------------
[1130]100         numror = 0                              ! define numror = 0 -> no restart file to read
[3]101         neuler = 0                              ! Set time-step indicator at nit000 (euler forward)
[1130]102         CALL day_init                           ! model calendar (using both namelist and restart infos)
[508]103         !                                       ! Initialization of ocean to zero
[3294]104         !   before fields      !       now fields     
105         sshb (:,:)   = 0._wp   ;   sshn (:,:)   = 0._wp
106         ub   (:,:,:) = 0._wp   ;   un   (:,:,:) = 0._wp
107         vb   (:,:,:) = 0._wp   ;   vn   (:,:,:) = 0._wp 
108         rotb (:,:,:) = 0._wp   ;   rotn (:,:,:) = 0._wp
109         hdivb(:,:,:) = 0._wp   ;   hdivn(:,:,:) = 0._wp
[508]110         !
[3]111         IF( cp_cfg == 'eel' ) THEN
[2104]112            CALL istate_eel                      ! EEL   configuration : start from pre-defined U,V T-S fields
[434]113         ELSEIF( cp_cfg == 'gyre' ) THEN         
[2104]114            CALL istate_gyre                     ! GYRE  configuration : start from pre-defined T-S fields
[4990]115        ELSEIF( cp_cfg == 'isomip' .OR. cp_cfg == 'isomip2') THEN
116            IF(lwp) WRITE(numout,*) 'Initialization of T+S for ISOMIP domain' 
117            tsn(:,:,:,jp_tem)=-1.9*tmask(:,:,:)          ! ISOMIP configuration : start from constant T+S fields
118            tsn(:,:,:,jp_sal)=34.4*tmask(:,:,:)
119            tsb(:,:,:,:)=tsn(:,:,:,:) 
[4245]120         ELSE                                    ! Initial T-S, U-V fields read in files
121            IF ( ln_tsd_init ) THEN              ! read 3D T and S data at nit000
122               CALL dta_tsd( nit000, tsb ) 
123               tsn(:,:,:,:) = tsb(:,:,:,:)
124               !
125            ELSE                                 ! Initial T-S fields defined analytically
126               CALL istate_t_s
127            ENDIF
128            IF ( ln_uvd_init .AND. lk_c1d ) THEN ! read 3D U and V data at nit000
129               CALL wrk_alloc( jpi, jpj, jpk, 2, zuvd )
130               CALL dta_uvd( nit000, zuvd )
131               ub(:,:,:) = zuvd(:,:,:,1) ;  un(:,:,:) = ub(:,:,:)
132               vb(:,:,:) = zuvd(:,:,:,2) ;  vn(:,:,:) = vb(:,:,:)
133               CALL wrk_dealloc( jpi, jpj, jpk, 2, zuvd )
134            ENDIF
[3]135         ENDIF
[2104]136         !
[4313]137         CALL eos( tsb, rhd, rhop, gdept_0(:,:,:) )        ! before potential and in situ densities
[2236]138#if ! defined key_c1d
[5120]139         IF( ln_zps .AND. .NOT. ln_isfcav)                                 &
140            &            CALL zps_hde    ( nit000, jpts, tsb, gtsu, gtsv,  &    ! Partial steps: before horizontal gradient
141            &                                            rhd, gru , grv    )  ! of t, s, rd at the last ocean level
142         IF( ln_zps .AND.       ln_isfcav)                                 &
[5189]143            &            CALL zps_hde_isf( nit000, jpts, tsb, gtsu, gtsv, gtui, gtvi,    &    ! Partial steps for top cell (ISF)
144            &                                            rhd, gru , grv , grui, grvi     )    ! of t, s, rd at the last ocean level
[2236]145#endif
[2148]146         !   
147         ! - ML - sshn could be modified by istate_eel, so that initialization of fse3t_b is done here
148         IF( lk_vvl ) THEN
149            DO jk = 1, jpk
150               fse3t_b(:,:,jk) = fse3t_n(:,:,jk)
151            ENDDO
152         ENDIF
153         !
[3]154      ENDIF
[1438]155      !
[2104]156      IF( lk_agrif ) THEN                  ! read free surface arrays in restart file
[593]157         IF( ln_rstart ) THEN
[3764]158            IF( lk_dynspg_flt )  THEN      ! read or initialize the following fields
159               !                           ! gcx, gcxb for agrif_opa_init
160               IF( sol_oce_alloc()  > 0 )   CALL ctl_stop('agrif sol_oce_alloc: allocation of arrays failed')
161               CALL flt_rst( nit000, 'READ' )
162            ENDIF
163         ENDIF                             ! explicit case not coded yet with AGRIF
[1200]164      ENDIF
[508]165      !
[4354]166      !
167      ! Initialize "now" and "before" barotropic velocities:
168      ! Do it whatever the free surface method, these arrays
169      ! being eventually used
170      !
171      !
172      un_b(:,:) = 0._wp ; vn_b(:,:) = 0._wp
173      ub_b(:,:) = 0._wp ; vb_b(:,:) = 0._wp
174      !
175      DO jk = 1, jpkm1
176         DO jj = 1, jpj
177            DO ji = 1, jpi
178               un_b(ji,jj) = un_b(ji,jj) + fse3u_n(ji,jj,jk) * un(ji,jj,jk) * umask(ji,jj,jk)
179               vn_b(ji,jj) = vn_b(ji,jj) + fse3v_n(ji,jj,jk) * vn(ji,jj,jk) * vmask(ji,jj,jk)
180               !
[4370]181               ub_b(ji,jj) = ub_b(ji,jj) + fse3u_b(ji,jj,jk) * ub(ji,jj,jk) * umask(ji,jj,jk)
182               vb_b(ji,jj) = vb_b(ji,jj) + fse3v_b(ji,jj,jk) * vb(ji,jj,jk) * vmask(ji,jj,jk)
[4354]183            END DO
184         END DO
185      END DO
186      !
[4370]187      un_b(:,:) = un_b(:,:) * hur  (:,:)
188      vn_b(:,:) = vn_b(:,:) * hvr  (:,:)
[4354]189      !
[4370]190      ub_b(:,:) = ub_b(:,:) * hur_b(:,:)
191      vb_b(:,:) = vb_b(:,:) * hvr_b(:,:)
[4354]192      !
193      !
[4990]194      IF( nn_timing == 1 )   CALL timing_stop('istate_init')
[3294]195      !
[3]196   END SUBROUTINE istate_init
197
[4990]198
[3294]199   SUBROUTINE istate_t_s
[3]200      !!---------------------------------------------------------------------
[3294]201      !!                  ***  ROUTINE istate_t_s  ***
[3]202      !!   
203      !! ** Purpose :   Intialization of the temperature field with an
204      !!      analytical profile or a file (i.e. in EEL configuration)
205      !!
[3294]206      !! ** Method  : - temperature: use Philander analytic profile
207      !!              - salinity   : use to a constant value 35.5
[3]208      !!
209      !! References :  Philander ???
210      !!----------------------------------------------------------------------
[3294]211      INTEGER  :: ji, jj, jk
212      REAL(wp) ::   zsal = 35.50
[3]213      !!----------------------------------------------------------------------
[508]214      !
[3]215      IF(lwp) WRITE(numout,*)
[3294]216      IF(lwp) WRITE(numout,*) 'istate_t_s : Philander s initial temperature profile'
217      IF(lwp) WRITE(numout,*) '~~~~~~~~~~   and constant salinity (',zsal,' psu)'
[2104]218      !
[3]219      DO jk = 1, jpk
[3294]220         tsn(:,:,jk,jp_tem) = (  ( ( 7.5 - 0. * ABS( gphit(:,:) )/30. ) * ( 1.-TANH((fsdept(:,:,jk)-80.)/30.) )   &
221            &                + 10. * ( 5000. - fsdept(:,:,jk) ) /5000.)  ) * tmask(:,:,jk)
222         tsb(:,:,jk,jp_tem) = tsn(:,:,jk,jp_tem)
[3]223      END DO
[3294]224      tsn(:,:,:,jp_sal) = zsal * tmask(:,:,:)
225      tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)
[2104]226      !
[3294]227   END SUBROUTINE istate_t_s
[3]228
[4990]229
[3]230   SUBROUTINE istate_eel
231      !!----------------------------------------------------------------------
232      !!                   ***  ROUTINE istate_eel  ***
233      !!
234      !! ** Purpose :   Initialization of the dynamics and tracers for EEL R5
235      !!      configuration (channel with or without a topographic bump)
236      !!
237      !! ** Method  : - set temprature field
238      !!              - set salinity field
239      !!              - set velocity field including horizontal divergence
240      !!                and relative vorticity fields
241      !!----------------------------------------------------------------------
242      USE divcur     ! hor. divergence & rel. vorticity      (div_cur routine)
[473]243      USE iom
[4990]244      !
[3]245      INTEGER  ::   inum              ! temporary logical unit
246      INTEGER  ::   ji, jj, jk        ! dummy loop indices
[479]247      INTEGER  ::   ijloc
[508]248      REAL(wp) ::   zh1, zh2, zslope, zcst, zfcor   ! temporary scalars
[2104]249      REAL(wp) ::   zt1  = 15._wp                   ! surface temperature value (EEL R5)
250      REAL(wp) ::   zt2  =  5._wp                   ! bottom  temperature value (EEL R5)
251      REAL(wp) ::   zsal = 35.0_wp                  ! constant salinity (EEL R2, R5 and R6)
252      REAL(wp) ::   zueel = 0.1_wp                  ! constant uniform zonal velocity (EEL R5)
[508]253      REAL(wp), DIMENSION(jpiglo,jpjglo) ::   zssh  ! initial ssh over the global domain
[3]254      !!----------------------------------------------------------------------
[4990]255      !
[3]256      SELECT CASE ( jp_cfg ) 
257         !                                              ! ====================
258         CASE ( 5 )                                     ! EEL R5 configuration
259            !                                           ! ====================
[2104]260            !
[3]261            ! set temperature field with a linear profile
262            ! -------------------------------------------
263            IF(lwp) WRITE(numout,*)
264            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: linear temperature profile'
265            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
[2104]266            !
[4292]267            zh1 = gdept_1d(  1  )
268            zh2 = gdept_1d(jpkm1)
[2104]269            !
[3]270            zslope = ( zt1 - zt2 ) / ( zh1 - zh2 )
271            zcst   = ( zt1 * ( zh1 - zh2) - ( zt1 - zt2 ) * zh1 ) / ( zh1 - zh2 )
[2104]272            !
[3]273            DO jk = 1, jpk
[3294]274               tsn(:,:,jk,jp_tem) = ( zt2 + zt1 * exp( - fsdept(:,:,jk) / 1000 ) ) * tmask(:,:,jk)
275               tsb(:,:,jk,jp_tem) = tsn(:,:,jk,jp_tem)
[3]276            END DO
[2104]277            !
[3294]278            IF(lwp) CALL prizre( tsn(:,:,:,jp_tem), jpi   , jpj   , jpk   , jpj/2 ,   &
279               &                             1     , jpi   , 5     , 1     , jpk   ,   &
280               &                             1     , 1.    , numout                  )
[2104]281            !
[3]282            ! set salinity field to a constant value
283            ! --------------------------------------
284            IF(lwp) WRITE(numout,*)
285            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal
286            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
[2104]287            !
[3294]288            tsn(:,:,:,jp_sal) = zsal * tmask(:,:,:)
289            tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)
[2104]290            !
[3]291            ! set the dynamics: U,V, hdiv, rot (and ssh if necessary)
292            ! ----------------
293            ! Start EEL5 configuration with barotropic geostrophic velocities
294            ! according the sshb and sshn SSH imposed.
[479]295            ! we assume a uniform grid (hence the use of e1t(1,1) for delta_y)
296            ! we use the Coriolis frequency at mid-channel.   
297            ub(:,:,:) = zueel * umask(:,:,:)
[3]298            un(:,:,:) = ub(:,:,:)
[479]299            ijloc = mj0(INT(jpjglo-1)/2)
300            zfcor = ff(1,ijloc)
[2104]301            !
[3]302            DO jj = 1, jpjglo
[479]303               zssh(:,jj) = - (FLOAT(jj)- FLOAT(jpjglo-1)/2.)*zueel*e1t(1,1)*zfcor/grav 
[3]304            END DO
[2104]305            !
[479]306            IF(lwp) THEN
307               WRITE(numout,*) ' Uniform zonal velocity for EEL R5:',zueel
308               WRITE(numout,*) ' Geostrophic SSH profile as a function of y:'
309               WRITE(numout,'(12(1x,f6.2))') zssh(1,:)
310            ENDIF
[2104]311            !
[3]312            DO jj = 1, nlcj
313               DO ji = 1, nlci
314                  sshb(ji,jj) = zssh( mig(ji) , mjg(jj) ) * tmask(ji,jj,1)
315               END DO
316            END DO
317            sshb(nlci+1:jpi,      :   ) = 0.e0      ! set to zero extra mpp columns
318            sshb(      :   ,nlcj+1:jpj) = 0.e0      ! set to zero extra mpp rows
[2104]319            !
[3]320            sshn(:,:) = sshb(:,:)                   ! set now ssh to the before value
[2104]321            !
[593]322            IF( nn_rstssh /= 0 ) THEN 
[2104]323               nn_rstssh = 0                        ! hand-made initilization of ssh
[593]324               CALL ctl_warn( 'istate_eel: force nn_rstssh = 0' )
[558]325            ENDIF
[2104]326            !
327            CALL div_cur( nit000 )                  ! horizontal divergence and relative vorticity (curl)
[3]328            ! N.B. the vertical velocity will be computed from the horizontal divergence field
329            ! in istate by a call to wzv routine
330
331
332            !                                     ! ==========================
333         CASE ( 2 , 6 )                           ! EEL R2 or R6 configuration
334            !                                     ! ==========================
[2104]335            !
[3]336            ! set temperature field with a NetCDF file
337            ! ----------------------------------------
338            IF(lwp) WRITE(numout,*)
339            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R2 or R6: read initial temperature in a NetCDF file'
340            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
[2104]341            !
[473]342            CALL iom_open ( 'eel.initemp', inum )
[3294]343            CALL iom_get ( inum, jpdom_data, 'initemp', tsb(:,:,:,jp_tem) ) ! read before temprature (tb)
[473]344            CALL iom_close( inum )
[2104]345            !
[3294]346            tsn(:,:,:,jp_tem) = tsb(:,:,:,jp_tem)                            ! set nox temperature to tb
[2104]347            !
[3294]348            IF(lwp) CALL prizre( tsn(:,:,:,jp_tem), jpi   , jpj   , jpk   , jpj/2 ,   &
349               &                            1     , jpi   , 5     , 1     , jpk   ,   &
350               &                            1     , 1.    , numout                  )
[2104]351            !
[3]352            ! set salinity field to a constant value
353            ! --------------------------------------
354            IF(lwp) WRITE(numout,*)
355            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal
356            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
[2104]357            !
[3294]358            tsn(:,:,:,jp_sal) = zsal * tmask(:,:,:)
359            tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)
[2104]360            !
[3]361            !                                    ! ===========================
362         CASE DEFAULT                            ! NONE existing configuration
363            !                                    ! ===========================
[473]364            WRITE(ctmp1,*) 'EEL with a ', jp_cfg,' km resolution is not coded'
365            CALL ctl_stop( ctmp1 )
[2104]366            !
[3]367      END SELECT
[2104]368      !
[3]369   END SUBROUTINE istate_eel
370
371
[93]372   SUBROUTINE istate_gyre
373      !!----------------------------------------------------------------------
374      !!                   ***  ROUTINE istate_gyre  ***
375      !!
376      !! ** Purpose :   Initialization of the dynamics and tracers for GYRE
377      !!      configuration (double gyre with rotated domain)
378      !!
379      !! ** Method  : - set temprature field
380      !!              - set salinity field
381      !!----------------------------------------------------------------------
[473]382      INTEGER :: ji, jj, jk  ! dummy loop indices
[508]383      INTEGER            ::   inum          ! temporary logical unit
384      INTEGER, PARAMETER ::   ntsinit = 0   ! (0/1) (analytical/input data files) T&S initialization
[93]385      !!----------------------------------------------------------------------
[4990]386      !
[434]387      SELECT CASE ( ntsinit)
[4990]388      !
[434]389      CASE ( 0 )                  ! analytical T/S profil deduced from LEVITUS
390         IF(lwp) WRITE(numout,*)
391         IF(lwp) WRITE(numout,*) 'istate_gyre : initial analytical T and S profil deduced from LEVITUS '
392         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
[4990]393         !
[434]394         DO jk = 1, jpk
395            DO jj = 1, jpj
396               DO ji = 1, jpi
[3294]397                  tsn(ji,jj,jk,jp_tem) = (  16. - 12. * TANH( (fsdept(ji,jj,jk) - 400) / 700 )         )   &
[434]398                       &           * (-TANH( (500-fsdept(ji,jj,jk)) / 150 ) + 1) / 2               &
399                       &       + (      15. * ( 1. - TANH( (fsdept(ji,jj,jk)-50.) / 1500.) )       &
400                       &                - 1.4 * TANH((fsdept(ji,jj,jk)-100.) / 100.)               &   
401                       &                + 7.  * (1500. - fsdept(ji,jj,jk)) / 1500.             )   & 
402                       &           * (-TANH( (fsdept(ji,jj,jk) - 500) / 150) + 1) / 2
[3294]403                  tsn(ji,jj,jk,jp_tem) = tsn(ji,jj,jk,jp_tem) * tmask(ji,jj,jk)
404                  tsb(ji,jj,jk,jp_tem) = tsn(ji,jj,jk,jp_tem)
[434]405
[3294]406                  tsn(ji,jj,jk,jp_sal) =  (  36.25 - 1.13 * TANH( (fsdept(ji,jj,jk) - 305) / 460 )  )  &
[434]407                     &              * (-TANH((500 - fsdept(ji,jj,jk)) / 150) + 1) / 2          &
408                     &          + (  35.55 + 1.25 * (5000. - fsdept(ji,jj,jk)) / 5000.         &
409                     &                - 1.62 * TANH( (fsdept(ji,jj,jk) - 60.  ) / 650. )       &
410                     &                + 0.2  * TANH( (fsdept(ji,jj,jk) - 35.  ) / 100. )       &
411                     &                + 0.2  * TANH( (fsdept(ji,jj,jk) - 1000.) / 5000.)    )  &
412                     &              * (-TANH((fsdept(ji,jj,jk) - 500) / 150) + 1) / 2 
[3294]413                  tsn(ji,jj,jk,jp_sal) = tsn(ji,jj,jk,jp_sal) * tmask(ji,jj,jk)
414                  tsb(ji,jj,jk,jp_sal) = tsn(ji,jj,jk,jp_sal)
[434]415               END DO
[93]416            END DO
417         END DO
[4990]418         !
[434]419      CASE ( 1 )                  ! T/S data fields read in dta_tem.nc/data_sal.nc files
420         IF(lwp) WRITE(numout,*)
421         IF(lwp) WRITE(numout,*) 'istate_gyre : initial T and S read from dta_tem.nc/data_sal.nc files'
422         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
423         IF(lwp) WRITE(numout,*) '              NetCDF FORMAT'
424
425         ! Read temperature field
426         ! ----------------------
[473]427         CALL iom_open ( 'data_tem', inum )
[3294]428         CALL iom_get ( inum, jpdom_data, 'votemper', tsn(:,:,:,jp_tem) ) 
[473]429         CALL iom_close( inum )
[434]430
[3294]431         tsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem) * tmask(:,:,:) 
432         tsb(:,:,:,jp_tem) = tsn(:,:,:,jp_tem)
[434]433
434         ! Read salinity field
435         ! -------------------
[473]436         CALL iom_open ( 'data_sal', inum )
[3294]437         CALL iom_get ( inum, jpdom_data, 'vosaline', tsn(:,:,:,jp_sal) ) 
[473]438         CALL iom_close( inum )
[434]439
[3294]440         tsn(:,:,:,jp_sal) = tsn(:,:,:,jp_sal) * tmask(:,:,:) 
441         tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)
[4990]442         !
[434]443      END SELECT
[4990]444      !
[93]445      IF(lwp) THEN
446         WRITE(numout,*)
447         WRITE(numout,*) '              Initial temperature and salinity profiles:'
[4292]448         WRITE(numout, "(9x,' level   gdept_1d   temperature   salinity   ')" )
449         WRITE(numout, "(10x, i4, 3f10.2)" ) ( jk, gdept_1d(jk), tsn(2,2,jk,jp_tem), tsn(2,2,jk,jp_sal), jk = 1, jpk )
[93]450      ENDIF
[4990]451      !
[93]452   END SUBROUTINE istate_gyre
453
[4990]454
[3]455   SUBROUTINE istate_uvg
456      !!----------------------------------------------------------------------
457      !!                  ***  ROUTINE istate_uvg  ***
458      !!
459      !! ** Purpose :   Compute the geostrophic velocities from (tn,sn) fields
460      !!
461      !! ** Method  :   Using the hydrostatic hypothesis the now hydrostatic
462      !!      pressure is computed by integrating the in-situ density from the
463      !!      surface to the bottom.
464      !!                 p=integral [ rau*g dz ]
465      !!----------------------------------------------------------------------
[359]466      USE dynspg          ! surface pressure gradient             (dyn_spg routine)
[3]467      USE divcur          ! hor. divergence & rel. vorticity      (div_cur routine)
468      USE lbclnk          ! ocean lateral boundary condition (or mpp link)
[4990]469      !
[3]470      INTEGER ::   ji, jj, jk        ! dummy loop indices
471      INTEGER ::   indic             ! ???
[508]472      REAL(wp) ::   zmsv, zphv, zmsu, zphu, zalfg     ! temporary scalars
[3294]473      REAL(wp), POINTER, DIMENSION(:,:,:) :: zprn
[3]474      !!----------------------------------------------------------------------
[3294]475      !
476      CALL wrk_alloc( jpi, jpj, jpk, zprn)
477      !
[3]478      IF(lwp) WRITE(numout,*) 
479      IF(lwp) WRITE(numout,*) 'istate_uvg : Start from Geostrophy'
480      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
481
482      ! Compute the now hydrostatic pressure
483      ! ------------------------------------
484
[15]485      zalfg = 0.5 * grav * rau0
[508]486     
487      zprn(:,:,1) = zalfg * fse3w(:,:,1) * ( 1 + rhd(:,:,1) )       ! Surface value
[3]488
[508]489      DO jk = 2, jpkm1                                              ! Vertical integration from the surface
[3]490         zprn(:,:,jk) = zprn(:,:,jk-1)   &
[359]491            &         + zalfg * fse3w(:,:,jk) * ( 2. + rhd(:,:,jk) + rhd(:,:,jk-1) )
[3]492      END DO 
493
494      ! Compute geostrophic balance
495      ! ---------------------------
496      DO jk = 1, jpkm1
497         DO jj = 2, jpjm1
498            DO ji = fs_2, fs_jpim1   ! vertor opt.
499               zmsv = 1. / MAX(  umask(ji-1,jj+1,jk) + umask(ji  ,jj+1,jk)   &
500                               + umask(ji-1,jj  ,jk) + umask(ji  ,jj  ,jk) , 1.  )
501               zphv = ( zprn(ji  ,jj+1,jk) - zprn(ji-1,jj+1,jk) ) * umask(ji-1,jj+1,jk) / e1u(ji-1,jj+1)   &
502                    + ( zprn(ji+1,jj+1,jk) - zprn(ji  ,jj+1,jk) ) * umask(ji  ,jj+1,jk) / e1u(ji  ,jj+1)   &
503                    + ( zprn(ji  ,jj  ,jk) - zprn(ji-1,jj  ,jk) ) * umask(ji-1,jj  ,jk) / e1u(ji-1,jj  )   &
504                    + ( zprn(ji+1,jj  ,jk) - zprn(ji  ,jj  ,jk) ) * umask(ji  ,jj  ,jk) / e1u(ji  ,jj  )
505               zphv = 1. / rau0 * zphv * zmsv * vmask(ji,jj,jk)
506
507               zmsu = 1. / MAX(  vmask(ji+1,jj  ,jk) + vmask(ji  ,jj  ,jk)   &
508                               + vmask(ji+1,jj-1,jk) + vmask(ji  ,jj-1,jk) , 1.  )
509               zphu = ( zprn(ji+1,jj+1,jk) - zprn(ji+1,jj  ,jk) ) * vmask(ji+1,jj  ,jk) / e2v(ji+1,jj  )   &
510                    + ( zprn(ji  ,jj+1,jk) - zprn(ji  ,jj  ,jk) ) * vmask(ji  ,jj  ,jk) / e2v(ji  ,jj  )   &
511                    + ( zprn(ji+1,jj  ,jk) - zprn(ji+1,jj-1,jk) ) * vmask(ji+1,jj-1,jk) / e2v(ji+1,jj-1)   &
512                    + ( zprn(ji  ,jj  ,jk) - zprn(ji  ,jj-1,jk) ) * vmask(ji  ,jj-1,jk) / e2v(ji  ,jj-1)
513               zphu = 1. / rau0 * zphu * zmsu * umask(ji,jj,jk)
514
515               ! Compute the geostrophic velocities
516               un(ji,jj,jk) = -2. * zphu / ( ff(ji,jj) + ff(ji  ,jj-1) )
517               vn(ji,jj,jk) =  2. * zphv / ( ff(ji,jj) + ff(ji-1,jj  ) )
518            END DO
519         END DO
520      END DO
521
522      IF(lwp) WRITE(numout,*) '         we force to zero bottom velocity'
523
524      ! Susbtract the bottom velocity (level jpk-1 for flat bottom case)
525      ! to have a zero bottom velocity
526
527      DO jk = 1, jpkm1
528         un(:,:,jk) = ( un(:,:,jk) - un(:,:,jpkm1) ) * umask(:,:,jk)
529         vn(:,:,jk) = ( vn(:,:,jk) - vn(:,:,jpkm1) ) * vmask(:,:,jk)
530      END DO
531
532      CALL lbc_lnk( un, 'U', -1. )
533      CALL lbc_lnk( vn, 'V', -1. )
534     
535      ub(:,:,:) = un(:,:,:)
536      vb(:,:,:) = vn(:,:,:)
537     
538      ! WARNING !!!!!
539      ! after initializing u and v, we need to calculate the initial streamfunction bsf.
540      ! Otherwise, only the trend will be computed and the model will blow up (inconsistency).
541      ! to do that, we call dyn_spg with a special trick:
[508]542      ! we fill ua and va with the velocities divided by dt, and the streamfunction will be brought to the
543      ! right value assuming the velocities have been set up in one time step.
544      ! we then set bsfd to zero (first guess for next step is d(psi)/dt = 0.)
545      !  sets up s false trend to calculate the barotropic streamfunction.
[3]546
547      ua(:,:,:) = ub(:,:,:) / rdt
548      va(:,:,:) = vb(:,:,:) / rdt
549
[359]550      ! calls dyn_spg. we assume euler time step, starting from rest.
[3]551      indic = 0
[359]552      CALL dyn_spg( nit000, indic )       ! surface pressure gradient
[3]553
554      ! the new velocity is ua*rdt
555
556      CALL lbc_lnk( ua, 'U', -1. )
557      CALL lbc_lnk( va, 'V', -1. )
558
559      ub(:,:,:) = ua(:,:,:) * rdt
560      vb(:,:,:) = va(:,:,:) * rdt
561      ua(:,:,:) = 0.e0
562      va(:,:,:) = 0.e0
563      un(:,:,:) = ub(:,:,:)
564      vn(:,:,:) = vb(:,:,:)
565       
566      ! Compute the divergence and curl
567
568      CALL div_cur( nit000 )            ! now horizontal divergence and curl
569
570      hdivb(:,:,:) = hdivn(:,:,:)       ! set the before to the now value
571      rotb (:,:,:) = rotn (:,:,:)       ! set the before to the now value
[508]572      !
[3294]573      CALL wrk_dealloc( jpi, jpj, jpk, zprn)
[2715]574      !
[3]575   END SUBROUTINE istate_uvg
576
577   !!=====================================================================
578END MODULE istate
Note: See TracBrowser for help on using the repository browser.