New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
istate.F90 in trunk/NEMOGCM/NEMO/OPA_SRC/DOM – NEMO

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

Last change on this file since 4528 was 4370, checked in by jchanut, 10 years ago

Time split cleaning / AMM12 restart / BDY and Agrif. See tickets #1228, #1227, #1225 and #1133

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