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/NEMOGCM/NEMO/OPA_SRC/DOM – NEMO

source: trunk/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 3693

Last change on this file since 3693 was 3693, checked in by rblod, 11 years ago

Fix Agrif and solver, see tocket #1030

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