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

source: branches/2013/dev_LOCEAN_CMCC_INGV_MERC_UKMO_2013/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 4245

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

dev_locean_cmcc_ingv_ukmo_merc : merge in the MERC_UKMO dev branch with trunk rev 4119

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