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

source: branches/2014/dev_r4621_NOC4_BDY_VERT_INTERP/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 5901

Last change on this file since 5901 was 5901, checked in by jamesharle, 8 years ago

merging branch with head of the trunk

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