New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
istate.F90 in branches/2014/dev_r4627_COMODO_UPW2D/NEMOGCM/CONFIG/UPW2D_eo_NOadv_vvl_EXP/MY_SRC – NEMO

source: branches/2014/dev_r4627_COMODO_UPW2D/NEMOGCM/CONFIG/UPW2D_eo_NOadv_vvl_EXP/MY_SRC/istate.F90 @ 4648

Last change on this file since 4648 was 4648, checked in by flavoni, 10 years ago

add new experience for cas test Upwelling, for WP item CNRS-7

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