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

source: branches/2013/dev_MERGE_2013/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 4292

Last change on this file since 4292 was 4292, checked in by cetlod, 10 years ago

dev_MERGE_2013 : 1st step of the merge, see ticket #1185

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