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

source: branches/2014/dev_CNRS0_NOC1_LDF/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 4616

Last change on this file since 4616 was 4616, checked in by gm, 10 years ago

#1260 : see the associated wiki page for explanation

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