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/NEMO/OPA_SRC – NEMO

source: trunk/NEMO/OPA_SRC/istate.F90 @ 1566

Last change on this file since 1566 was 1566, checked in by rblod, 15 years ago

Cosmetic changes: suppress useless variables and code review of the code changed when suppressing rigid-lid, see ticket #508

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