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

Last change on this file since 591 was 558, checked in by opalod, 18 years ago

nemo_v1_bugfix_070: SM+RB+GM: add nn_rstssh for hand made initilization of ssh or not (1/0)

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