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

source: branches/2013/dev_r3858_NOC_ZTC/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 3862

Last change on this file since 3862 was 3862, checked in by acc, 11 years ago

Branch 2013/dev_r3858_NOC_ZTC #863. Change names of gdept_0, gdepw_0, e3t_0, e3w_0 to gdept_1d, gdepw_1d, e3t_1d, e3w_1d respectively in preparation for the introduction of time invariant 3d arrays using the _0 suffix

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