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 @ 5332

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

ISOMIP config code removed (mv to unsupported configuration); correction of bug about initialisation of top friction (ISF only), about coastline modification if ISF activated and top and bottom e3uw if ISF.

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