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_r5094_UKMO_ISFCLEAN/NEMOGCM/NEMO/OPA_SRC/DOM – NEMO

source: branches/2015/dev_r5094_UKMO_ISFCLEAN/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 5111

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

add some missing if ln_isfcav, test of option compatibility with ln_isfcav + small documentation

  • Property svn:keywords set to Id
File size: 27.5 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
[5104]139         IF( ln_zps .AND. .NOT. ln_isfcav)                                 &
140            &            CALL zps_hde    ( nit000, jpts, tsb, gtsu, gtsv,  &    ! Partial steps: before horizontal gradient
[5098]141            &                                            rhd, gru , grv    )  ! of t, s, rd at the last ocean level
[5111]142         IF( ln_zps .AND.       ln_isfcav)                                 &
[5098]143            &            CALL zps_hde_isf( nit000, jpts, tsb, gtsu, gtsv,  &    ! Partial steps for top cell (ISF)
144            &                                            rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv ,   &
145            &                                     gtui, gtvi, grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi    ) ! of t, s, rd at the last ocean level
[2236]146#endif
[2148]147         !   
148         ! - ML - sshn could be modified by istate_eel, so that initialization of fse3t_b is done here
149         IF( lk_vvl ) THEN
150            DO jk = 1, jpk
151               fse3t_b(:,:,jk) = fse3t_n(:,:,jk)
152            ENDDO
153         ENDIF
154         !
[3]155      ENDIF
[1438]156      !
[2104]157      IF( lk_agrif ) THEN                  ! read free surface arrays in restart file
[593]158         IF( ln_rstart ) THEN
[3764]159            IF( lk_dynspg_flt )  THEN      ! read or initialize the following fields
160               !                           ! gcx, gcxb for agrif_opa_init
161               IF( sol_oce_alloc()  > 0 )   CALL ctl_stop('agrif sol_oce_alloc: allocation of arrays failed')
162               CALL flt_rst( nit000, 'READ' )
163            ENDIF
164         ENDIF                             ! explicit case not coded yet with AGRIF
[1200]165      ENDIF
[508]166      !
[4354]167      !
168      ! Initialize "now" and "before" barotropic velocities:
169      ! Do it whatever the free surface method, these arrays
170      ! being eventually used
171      !
172      !
173      un_b(:,:) = 0._wp ; vn_b(:,:) = 0._wp
174      ub_b(:,:) = 0._wp ; vb_b(:,:) = 0._wp
175      !
176      DO jk = 1, jpkm1
177         DO jj = 1, jpj
178            DO ji = 1, jpi
179               un_b(ji,jj) = un_b(ji,jj) + fse3u_n(ji,jj,jk) * un(ji,jj,jk) * umask(ji,jj,jk)
180               vn_b(ji,jj) = vn_b(ji,jj) + fse3v_n(ji,jj,jk) * vn(ji,jj,jk) * vmask(ji,jj,jk)
181               !
[4370]182               ub_b(ji,jj) = ub_b(ji,jj) + fse3u_b(ji,jj,jk) * ub(ji,jj,jk) * umask(ji,jj,jk)
183               vb_b(ji,jj) = vb_b(ji,jj) + fse3v_b(ji,jj,jk) * vb(ji,jj,jk) * vmask(ji,jj,jk)
[4354]184            END DO
185         END DO
186      END DO
187      !
[4370]188      un_b(:,:) = un_b(:,:) * hur  (:,:)
189      vn_b(:,:) = vn_b(:,:) * hvr  (:,:)
[4354]190      !
[4370]191      ub_b(:,:) = ub_b(:,:) * hur_b(:,:)
192      vb_b(:,:) = vb_b(:,:) * hvr_b(:,:)
[4354]193      !
194      !
[4990]195      IF( nn_timing == 1 )   CALL timing_stop('istate_init')
[3294]196      !
[3]197   END SUBROUTINE istate_init
198
[4990]199
[3294]200   SUBROUTINE istate_t_s
[3]201      !!---------------------------------------------------------------------
[3294]202      !!                  ***  ROUTINE istate_t_s  ***
[3]203      !!   
204      !! ** Purpose :   Intialization of the temperature field with an
205      !!      analytical profile or a file (i.e. in EEL configuration)
206      !!
[3294]207      !! ** Method  : - temperature: use Philander analytic profile
208      !!              - salinity   : use to a constant value 35.5
[3]209      !!
210      !! References :  Philander ???
211      !!----------------------------------------------------------------------
[3294]212      INTEGER  :: ji, jj, jk
213      REAL(wp) ::   zsal = 35.50
[3]214      !!----------------------------------------------------------------------
[508]215      !
[3]216      IF(lwp) WRITE(numout,*)
[3294]217      IF(lwp) WRITE(numout,*) 'istate_t_s : Philander s initial temperature profile'
218      IF(lwp) WRITE(numout,*) '~~~~~~~~~~   and constant salinity (',zsal,' psu)'
[2104]219      !
[3]220      DO jk = 1, jpk
[3294]221         tsn(:,:,jk,jp_tem) = (  ( ( 7.5 - 0. * ABS( gphit(:,:) )/30. ) * ( 1.-TANH((fsdept(:,:,jk)-80.)/30.) )   &
222            &                + 10. * ( 5000. - fsdept(:,:,jk) ) /5000.)  ) * tmask(:,:,jk)
223         tsb(:,:,jk,jp_tem) = tsn(:,:,jk,jp_tem)
[3]224      END DO
[3294]225      tsn(:,:,:,jp_sal) = zsal * tmask(:,:,:)
226      tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)
[2104]227      !
[3294]228   END SUBROUTINE istate_t_s
[3]229
[4990]230
[3]231   SUBROUTINE istate_eel
232      !!----------------------------------------------------------------------
233      !!                   ***  ROUTINE istate_eel  ***
234      !!
235      !! ** Purpose :   Initialization of the dynamics and tracers for EEL R5
236      !!      configuration (channel with or without a topographic bump)
237      !!
238      !! ** Method  : - set temprature field
239      !!              - set salinity field
240      !!              - set velocity field including horizontal divergence
241      !!                and relative vorticity fields
242      !!----------------------------------------------------------------------
243      USE divcur     ! hor. divergence & rel. vorticity      (div_cur routine)
[473]244      USE iom
[4990]245      !
[3]246      INTEGER  ::   inum              ! temporary logical unit
247      INTEGER  ::   ji, jj, jk        ! dummy loop indices
[479]248      INTEGER  ::   ijloc
[508]249      REAL(wp) ::   zh1, zh2, zslope, zcst, zfcor   ! temporary scalars
[2104]250      REAL(wp) ::   zt1  = 15._wp                   ! surface temperature value (EEL R5)
251      REAL(wp) ::   zt2  =  5._wp                   ! bottom  temperature value (EEL R5)
252      REAL(wp) ::   zsal = 35.0_wp                  ! constant salinity (EEL R2, R5 and R6)
253      REAL(wp) ::   zueel = 0.1_wp                  ! constant uniform zonal velocity (EEL R5)
[508]254      REAL(wp), DIMENSION(jpiglo,jpjglo) ::   zssh  ! initial ssh over the global domain
[3]255      !!----------------------------------------------------------------------
[4990]256      !
[3]257      SELECT CASE ( jp_cfg ) 
258         !                                              ! ====================
259         CASE ( 5 )                                     ! EEL R5 configuration
260            !                                           ! ====================
[2104]261            !
[3]262            ! set temperature field with a linear profile
263            ! -------------------------------------------
264            IF(lwp) WRITE(numout,*)
265            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: linear temperature profile'
266            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
[2104]267            !
[4292]268            zh1 = gdept_1d(  1  )
269            zh2 = gdept_1d(jpkm1)
[2104]270            !
[3]271            zslope = ( zt1 - zt2 ) / ( zh1 - zh2 )
272            zcst   = ( zt1 * ( zh1 - zh2) - ( zt1 - zt2 ) * zh1 ) / ( zh1 - zh2 )
[2104]273            !
[3]274            DO jk = 1, jpk
[3294]275               tsn(:,:,jk,jp_tem) = ( zt2 + zt1 * exp( - fsdept(:,:,jk) / 1000 ) ) * tmask(:,:,jk)
276               tsb(:,:,jk,jp_tem) = tsn(:,:,jk,jp_tem)
[3]277            END DO
[2104]278            !
[3294]279            IF(lwp) CALL prizre( tsn(:,:,:,jp_tem), jpi   , jpj   , jpk   , jpj/2 ,   &
280               &                             1     , jpi   , 5     , 1     , jpk   ,   &
281               &                             1     , 1.    , numout                  )
[2104]282            !
[3]283            ! set salinity field to a constant value
284            ! --------------------------------------
285            IF(lwp) WRITE(numout,*)
286            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal
287            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
[2104]288            !
[3294]289            tsn(:,:,:,jp_sal) = zsal * tmask(:,:,:)
290            tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)
[2104]291            !
[3]292            ! set the dynamics: U,V, hdiv, rot (and ssh if necessary)
293            ! ----------------
294            ! Start EEL5 configuration with barotropic geostrophic velocities
295            ! according the sshb and sshn SSH imposed.
[479]296            ! we assume a uniform grid (hence the use of e1t(1,1) for delta_y)
297            ! we use the Coriolis frequency at mid-channel.   
298            ub(:,:,:) = zueel * umask(:,:,:)
[3]299            un(:,:,:) = ub(:,:,:)
[479]300            ijloc = mj0(INT(jpjglo-1)/2)
301            zfcor = ff(1,ijloc)
[2104]302            !
[3]303            DO jj = 1, jpjglo
[479]304               zssh(:,jj) = - (FLOAT(jj)- FLOAT(jpjglo-1)/2.)*zueel*e1t(1,1)*zfcor/grav 
[3]305            END DO
[2104]306            !
[479]307            IF(lwp) THEN
308               WRITE(numout,*) ' Uniform zonal velocity for EEL R5:',zueel
309               WRITE(numout,*) ' Geostrophic SSH profile as a function of y:'
310               WRITE(numout,'(12(1x,f6.2))') zssh(1,:)
311            ENDIF
[2104]312            !
[3]313            DO jj = 1, nlcj
314               DO ji = 1, nlci
315                  sshb(ji,jj) = zssh( mig(ji) , mjg(jj) ) * tmask(ji,jj,1)
316               END DO
317            END DO
318            sshb(nlci+1:jpi,      :   ) = 0.e0      ! set to zero extra mpp columns
319            sshb(      :   ,nlcj+1:jpj) = 0.e0      ! set to zero extra mpp rows
[2104]320            !
[3]321            sshn(:,:) = sshb(:,:)                   ! set now ssh to the before value
[2104]322            !
[593]323            IF( nn_rstssh /= 0 ) THEN 
[2104]324               nn_rstssh = 0                        ! hand-made initilization of ssh
[593]325               CALL ctl_warn( 'istate_eel: force nn_rstssh = 0' )
[558]326            ENDIF
[2104]327            !
328            CALL div_cur( nit000 )                  ! horizontal divergence and relative vorticity (curl)
[3]329            ! N.B. the vertical velocity will be computed from the horizontal divergence field
330            ! in istate by a call to wzv routine
331
332
333            !                                     ! ==========================
334         CASE ( 2 , 6 )                           ! EEL R2 or R6 configuration
335            !                                     ! ==========================
[2104]336            !
[3]337            ! set temperature field with a NetCDF file
338            ! ----------------------------------------
339            IF(lwp) WRITE(numout,*)
340            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R2 or R6: read initial temperature in a NetCDF file'
341            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
[2104]342            !
[473]343            CALL iom_open ( 'eel.initemp', inum )
[3294]344            CALL iom_get ( inum, jpdom_data, 'initemp', tsb(:,:,:,jp_tem) ) ! read before temprature (tb)
[473]345            CALL iom_close( inum )
[2104]346            !
[3294]347            tsn(:,:,:,jp_tem) = tsb(:,:,:,jp_tem)                            ! set nox temperature to tb
[2104]348            !
[3294]349            IF(lwp) CALL prizre( tsn(:,:,:,jp_tem), jpi   , jpj   , jpk   , jpj/2 ,   &
350               &                            1     , jpi   , 5     , 1     , jpk   ,   &
351               &                            1     , 1.    , numout                  )
[2104]352            !
[3]353            ! set salinity field to a constant value
354            ! --------------------------------------
355            IF(lwp) WRITE(numout,*)
356            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal
357            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
[2104]358            !
[3294]359            tsn(:,:,:,jp_sal) = zsal * tmask(:,:,:)
360            tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)
[2104]361            !
[3]362            !                                    ! ===========================
363         CASE DEFAULT                            ! NONE existing configuration
364            !                                    ! ===========================
[473]365            WRITE(ctmp1,*) 'EEL with a ', jp_cfg,' km resolution is not coded'
366            CALL ctl_stop( ctmp1 )
[2104]367            !
[3]368      END SELECT
[2104]369      !
[3]370   END SUBROUTINE istate_eel
371
372
[93]373   SUBROUTINE istate_gyre
374      !!----------------------------------------------------------------------
375      !!                   ***  ROUTINE istate_gyre  ***
376      !!
377      !! ** Purpose :   Initialization of the dynamics and tracers for GYRE
378      !!      configuration (double gyre with rotated domain)
379      !!
380      !! ** Method  : - set temprature field
381      !!              - set salinity field
382      !!----------------------------------------------------------------------
[473]383      INTEGER :: ji, jj, jk  ! dummy loop indices
[508]384      INTEGER            ::   inum          ! temporary logical unit
385      INTEGER, PARAMETER ::   ntsinit = 0   ! (0/1) (analytical/input data files) T&S initialization
[93]386      !!----------------------------------------------------------------------
[4990]387      !
[434]388      SELECT CASE ( ntsinit)
[4990]389      !
[434]390      CASE ( 0 )                  ! analytical T/S profil deduced from LEVITUS
391         IF(lwp) WRITE(numout,*)
392         IF(lwp) WRITE(numout,*) 'istate_gyre : initial analytical T and S profil deduced from LEVITUS '
393         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
[4990]394         !
[434]395         DO jk = 1, jpk
396            DO jj = 1, jpj
397               DO ji = 1, jpi
[3294]398                  tsn(ji,jj,jk,jp_tem) = (  16. - 12. * TANH( (fsdept(ji,jj,jk) - 400) / 700 )         )   &
[434]399                       &           * (-TANH( (500-fsdept(ji,jj,jk)) / 150 ) + 1) / 2               &
400                       &       + (      15. * ( 1. - TANH( (fsdept(ji,jj,jk)-50.) / 1500.) )       &
401                       &                - 1.4 * TANH((fsdept(ji,jj,jk)-100.) / 100.)               &   
402                       &                + 7.  * (1500. - fsdept(ji,jj,jk)) / 1500.             )   & 
403                       &           * (-TANH( (fsdept(ji,jj,jk) - 500) / 150) + 1) / 2
[3294]404                  tsn(ji,jj,jk,jp_tem) = tsn(ji,jj,jk,jp_tem) * tmask(ji,jj,jk)
405                  tsb(ji,jj,jk,jp_tem) = tsn(ji,jj,jk,jp_tem)
[434]406
[3294]407                  tsn(ji,jj,jk,jp_sal) =  (  36.25 - 1.13 * TANH( (fsdept(ji,jj,jk) - 305) / 460 )  )  &
[434]408                     &              * (-TANH((500 - fsdept(ji,jj,jk)) / 150) + 1) / 2          &
409                     &          + (  35.55 + 1.25 * (5000. - fsdept(ji,jj,jk)) / 5000.         &
410                     &                - 1.62 * TANH( (fsdept(ji,jj,jk) - 60.  ) / 650. )       &
411                     &                + 0.2  * TANH( (fsdept(ji,jj,jk) - 35.  ) / 100. )       &
412                     &                + 0.2  * TANH( (fsdept(ji,jj,jk) - 1000.) / 5000.)    )  &
413                     &              * (-TANH((fsdept(ji,jj,jk) - 500) / 150) + 1) / 2 
[3294]414                  tsn(ji,jj,jk,jp_sal) = tsn(ji,jj,jk,jp_sal) * tmask(ji,jj,jk)
415                  tsb(ji,jj,jk,jp_sal) = tsn(ji,jj,jk,jp_sal)
[434]416               END DO
[93]417            END DO
418         END DO
[4990]419         !
[434]420      CASE ( 1 )                  ! T/S data fields read in dta_tem.nc/data_sal.nc files
421         IF(lwp) WRITE(numout,*)
422         IF(lwp) WRITE(numout,*) 'istate_gyre : initial T and S read from dta_tem.nc/data_sal.nc files'
423         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
424         IF(lwp) WRITE(numout,*) '              NetCDF FORMAT'
425
426         ! Read temperature field
427         ! ----------------------
[473]428         CALL iom_open ( 'data_tem', inum )
[3294]429         CALL iom_get ( inum, jpdom_data, 'votemper', tsn(:,:,:,jp_tem) ) 
[473]430         CALL iom_close( inum )
[434]431
[3294]432         tsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem) * tmask(:,:,:) 
433         tsb(:,:,:,jp_tem) = tsn(:,:,:,jp_tem)
[434]434
435         ! Read salinity field
436         ! -------------------
[473]437         CALL iom_open ( 'data_sal', inum )
[3294]438         CALL iom_get ( inum, jpdom_data, 'vosaline', tsn(:,:,:,jp_sal) ) 
[473]439         CALL iom_close( inum )
[434]440
[3294]441         tsn(:,:,:,jp_sal) = tsn(:,:,:,jp_sal) * tmask(:,:,:) 
442         tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)
[4990]443         !
[434]444      END SELECT
[4990]445      !
[93]446      IF(lwp) THEN
447         WRITE(numout,*)
448         WRITE(numout,*) '              Initial temperature and salinity profiles:'
[4292]449         WRITE(numout, "(9x,' level   gdept_1d   temperature   salinity   ')" )
450         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]451      ENDIF
[4990]452      !
[93]453   END SUBROUTINE istate_gyre
454
[4990]455
[3]456   SUBROUTINE istate_uvg
457      !!----------------------------------------------------------------------
458      !!                  ***  ROUTINE istate_uvg  ***
459      !!
460      !! ** Purpose :   Compute the geostrophic velocities from (tn,sn) fields
461      !!
462      !! ** Method  :   Using the hydrostatic hypothesis the now hydrostatic
463      !!      pressure is computed by integrating the in-situ density from the
464      !!      surface to the bottom.
465      !!                 p=integral [ rau*g dz ]
466      !!----------------------------------------------------------------------
[359]467      USE dynspg          ! surface pressure gradient             (dyn_spg routine)
[3]468      USE divcur          ! hor. divergence & rel. vorticity      (div_cur routine)
469      USE lbclnk          ! ocean lateral boundary condition (or mpp link)
[4990]470      !
[3]471      INTEGER ::   ji, jj, jk        ! dummy loop indices
472      INTEGER ::   indic             ! ???
[508]473      REAL(wp) ::   zmsv, zphv, zmsu, zphu, zalfg     ! temporary scalars
[3294]474      REAL(wp), POINTER, DIMENSION(:,:,:) :: zprn
[3]475      !!----------------------------------------------------------------------
[3294]476      !
477      CALL wrk_alloc( jpi, jpj, jpk, zprn)
478      !
[3]479      IF(lwp) WRITE(numout,*) 
480      IF(lwp) WRITE(numout,*) 'istate_uvg : Start from Geostrophy'
481      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
482
483      ! Compute the now hydrostatic pressure
484      ! ------------------------------------
485
[15]486      zalfg = 0.5 * grav * rau0
[508]487     
488      zprn(:,:,1) = zalfg * fse3w(:,:,1) * ( 1 + rhd(:,:,1) )       ! Surface value
[3]489
[508]490      DO jk = 2, jpkm1                                              ! Vertical integration from the surface
[3]491         zprn(:,:,jk) = zprn(:,:,jk-1)   &
[359]492            &         + zalfg * fse3w(:,:,jk) * ( 2. + rhd(:,:,jk) + rhd(:,:,jk-1) )
[3]493      END DO 
494
495      ! Compute geostrophic balance
496      ! ---------------------------
497      DO jk = 1, jpkm1
498         DO jj = 2, jpjm1
499            DO ji = fs_2, fs_jpim1   ! vertor opt.
500               zmsv = 1. / MAX(  umask(ji-1,jj+1,jk) + umask(ji  ,jj+1,jk)   &
501                               + umask(ji-1,jj  ,jk) + umask(ji  ,jj  ,jk) , 1.  )
502               zphv = ( zprn(ji  ,jj+1,jk) - zprn(ji-1,jj+1,jk) ) * umask(ji-1,jj+1,jk) / e1u(ji-1,jj+1)   &
503                    + ( zprn(ji+1,jj+1,jk) - zprn(ji  ,jj+1,jk) ) * umask(ji  ,jj+1,jk) / e1u(ji  ,jj+1)   &
504                    + ( zprn(ji  ,jj  ,jk) - zprn(ji-1,jj  ,jk) ) * umask(ji-1,jj  ,jk) / e1u(ji-1,jj  )   &
505                    + ( zprn(ji+1,jj  ,jk) - zprn(ji  ,jj  ,jk) ) * umask(ji  ,jj  ,jk) / e1u(ji  ,jj  )
506               zphv = 1. / rau0 * zphv * zmsv * vmask(ji,jj,jk)
507
508               zmsu = 1. / MAX(  vmask(ji+1,jj  ,jk) + vmask(ji  ,jj  ,jk)   &
509                               + vmask(ji+1,jj-1,jk) + vmask(ji  ,jj-1,jk) , 1.  )
510               zphu = ( zprn(ji+1,jj+1,jk) - zprn(ji+1,jj  ,jk) ) * vmask(ji+1,jj  ,jk) / e2v(ji+1,jj  )   &
511                    + ( zprn(ji  ,jj+1,jk) - zprn(ji  ,jj  ,jk) ) * vmask(ji  ,jj  ,jk) / e2v(ji  ,jj  )   &
512                    + ( zprn(ji+1,jj  ,jk) - zprn(ji+1,jj-1,jk) ) * vmask(ji+1,jj-1,jk) / e2v(ji+1,jj-1)   &
513                    + ( zprn(ji  ,jj  ,jk) - zprn(ji  ,jj-1,jk) ) * vmask(ji  ,jj-1,jk) / e2v(ji  ,jj-1)
514               zphu = 1. / rau0 * zphu * zmsu * umask(ji,jj,jk)
515
516               ! Compute the geostrophic velocities
517               un(ji,jj,jk) = -2. * zphu / ( ff(ji,jj) + ff(ji  ,jj-1) )
518               vn(ji,jj,jk) =  2. * zphv / ( ff(ji,jj) + ff(ji-1,jj  ) )
519            END DO
520         END DO
521      END DO
522
523      IF(lwp) WRITE(numout,*) '         we force to zero bottom velocity'
524
525      ! Susbtract the bottom velocity (level jpk-1 for flat bottom case)
526      ! to have a zero bottom velocity
527
528      DO jk = 1, jpkm1
529         un(:,:,jk) = ( un(:,:,jk) - un(:,:,jpkm1) ) * umask(:,:,jk)
530         vn(:,:,jk) = ( vn(:,:,jk) - vn(:,:,jpkm1) ) * vmask(:,:,jk)
531      END DO
532
533      CALL lbc_lnk( un, 'U', -1. )
534      CALL lbc_lnk( vn, 'V', -1. )
535     
536      ub(:,:,:) = un(:,:,:)
537      vb(:,:,:) = vn(:,:,:)
538     
539      ! WARNING !!!!!
540      ! after initializing u and v, we need to calculate the initial streamfunction bsf.
541      ! Otherwise, only the trend will be computed and the model will blow up (inconsistency).
542      ! to do that, we call dyn_spg with a special trick:
[508]543      ! we fill ua and va with the velocities divided by dt, and the streamfunction will be brought to the
544      ! right value assuming the velocities have been set up in one time step.
545      ! we then set bsfd to zero (first guess for next step is d(psi)/dt = 0.)
546      !  sets up s false trend to calculate the barotropic streamfunction.
[3]547
548      ua(:,:,:) = ub(:,:,:) / rdt
549      va(:,:,:) = vb(:,:,:) / rdt
550
[359]551      ! calls dyn_spg. we assume euler time step, starting from rest.
[3]552      indic = 0
[359]553      CALL dyn_spg( nit000, indic )       ! surface pressure gradient
[3]554
555      ! the new velocity is ua*rdt
556
557      CALL lbc_lnk( ua, 'U', -1. )
558      CALL lbc_lnk( va, 'V', -1. )
559
560      ub(:,:,:) = ua(:,:,:) * rdt
561      vb(:,:,:) = va(:,:,:) * rdt
562      ua(:,:,:) = 0.e0
563      va(:,:,:) = 0.e0
564      un(:,:,:) = ub(:,:,:)
565      vn(:,:,:) = vb(:,:,:)
566       
567      ! Compute the divergence and curl
568
569      CALL div_cur( nit000 )            ! now horizontal divergence and curl
570
571      hdivb(:,:,:) = hdivn(:,:,:)       ! set the before to the now value
572      rotb (:,:,:) = rotn (:,:,:)       ! set the before to the now value
[508]573      !
[3294]574      CALL wrk_dealloc( jpi, jpj, jpk, zprn)
[2715]575      !
[3]576   END SUBROUTINE istate_uvg
577
578   !!=====================================================================
579END MODULE istate
Note: See TracBrowser for help on using the repository browser.