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

Last change on this file since 703 was 699, checked in by smasson, 17 years ago

insert revision Id

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