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

source: branches/2011/dev_NEMO_MERGE_2011/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 3115

Last change on this file since 3115 was 2977, checked in by cetlod, 13 years ago

Add in branch 2011/dev_LOCEAN_2011 changes from 2011/dev_r2787_PISCES_improvment, 2011/dev_r2787_LOCEAN_offline_fldread and 2011/dev_r2787_LOCEAN3_TRA_TRP branches, see ticket #877

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