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

source: branches/2015/dev_r5721_CNRS9_NOC3_LDF/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 5770

Last change on this file since 5770 was 5770, checked in by gm, 9 years ago

#1593: LDF-ADV, step II.2: phasing the improvements/simplifications of advective tracer trend (see wiki page)

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