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

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

nemo_v2_update_001 : CT : - add non linear free surface (variable volume) with new cpp key key_vvl

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