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 @ 5766

Last change on this file since 5766 was 5766, checked in by cetlod, 9 years ago

LDF: phasing the improvements/simplifications of TOP component

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