New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
istate.F90 in branches/nemo_v3_3_beta/NEMOGCM/NEMO/OPA_SRC – NEMO

source: branches/nemo_v3_3_beta/NEMOGCM/NEMO/OPA_SRC/istate.F90 @ 2281

Last change on this file since 2281 was 2281, checked in by smasson, 14 years ago

set proper svn properties to all files...

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