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

Last change on this file since 508 was 508, checked in by opalod, 17 years ago

nemo_v1_update_071:RB: add iom for restart and reorganization of restart

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