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

source: branches/2012/dev_v3_4_STABLE_2012/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 5169

Last change on this file since 5169 was 5169, checked in by pabouttier, 9 years ago

Fix missing rn2b initialisation (direct model and TAM), described in Ticket #1498

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