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

source: trunk/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 5163

Last change on this file since 5163 was 5163, checked in by gm, 9 years ago

trunk: #1498 bug correction on rn2b initialization

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