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

source: branches/UKMO/dev_r5107_hadgem3_cplseq/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 5477

Last change on this file since 5477 was 5477, checked in by cguiavarch, 9 years ago

Clear svn keywords from UKMO/dev_r5107_hadgem3_cplseq

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