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

source: branches/2012/dev_v3_4_STABLE_2012/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 5168

Last change on this file since 5168 was 5168, checked in by pabouttier, 7 years ago

Fixed several bugs described in Ticket #1360

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