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

source: branches/2015/dev_r5151_UKMO_ISF/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 5621

Last change on this file since 5621 was 5621, checked in by mathiot, 9 years ago

UKMO_ISF: upgrade to NEMO_3.6_STABLE (r5554)

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