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

Last change on this file since 3562 was 3562, checked in by rblod, 11 years ago

Fix vvl and flux formulae, see ticket #952

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