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

source: branches/2013/dev_MERGE_2013/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 4354

Last change on this file since 4354 was 4354, checked in by jchanut, 10 years ago

Restore AGRIF and BDY compatibility, see ticket #1133

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