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

source: branches/UKMO/dev_r5518_medusa_fix_restart/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 7865

Last change on this file since 7865 was 7865, checked in by marc, 7 years ago

Adding extra arguments to ZPS_HDE

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