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

source: branches/2012/dev_v3_4_STABLE_2012/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 5196

Last change on this file since 5196 was 5196, checked in by pabouttier, 9 years ago

Add SEABASS reference configuration for this NEMO version for now; See Ticket #1505

  • Property svn:keywords set to Id
File size: 27.3 KB
Line 
1MODULE istate
2   !!======================================================================
3   !!                     ***  MODULE  istate  ***
4   !! Ocean state   :  initial state setting
5   !!=====================================================================
6   !! History :  OPA  !  1989-12  (P. Andrich)  Original code
7   !!            5.0  !  1991-11  (G. Madec)  rewritting
8   !!            6.0  !  1996-01  (G. Madec)  terrain following coordinates
9   !!            8.0  !  2001-09  (M. Levy, M. Ben Jelloul)  istate_eel
10   !!            8.0  !  2001-09  (M. Levy, M. Ben Jelloul)  istate_uvg
11   !!   NEMO     1.0  !  2003-08  (G. Madec, C. Talandier)  F90: Free form, modules + EEL R5
12   !!             -   !  2004-05  (A. Koch-Larrouy)  istate_gyre
13   !!            2.0  !  2006-07  (S. Masson)  distributed restart using iom
14   !!            3.3  !  2010-10  (C. Ethe) merge TRC-TRA
15   !!            3.4  !  2011-04  (G. Madec) Merge of dtatem and dtasal & suppression of tb,tn/sb,sn
16   !!----------------------------------------------------------------------
17
18   !!----------------------------------------------------------------------
19   !!   istate_init   : initial state setting
20   !!   istate_tem    : analytical profile for initial Temperature
21   !!   istate_sal    : analytical profile for initial Salinity
22   !!   istate_eel    : initial state setting of EEL R5 configuration
23   !!   istate_gyre   : initial state setting of GYRE configuration
24   !!   istate_uvg    : initial velocity in geostropic balance
25   !!----------------------------------------------------------------------
26   USE oce             ! ocean dynamics and active tracers
27   USE dom_oce         ! ocean space and time domain
28   USE daymod          ! calendar
29   USE eosbn2          ! eq. of state, Brunt Vaisala frequency (eos     routine)
30   USE ldftra_oce      ! ocean active tracers: lateral physics
31   USE zdf_oce         ! ocean vertical physics
32   USE phycst          ! physical constants
33   USE dtatsd          ! data temperature and salinity   (dta_tsd routine)
34   USE restart         ! ocean restart                   (rst_read routine)
35   USE in_out_manager  ! I/O manager
36   USE iom             ! I/O library
37   USE zpshde          ! partial step: hor. derivative (zps_hde routine)
38   USE eosbn2          ! equation of state            (eos bn2 routine)
39   USE domvvl          ! varying vertical mesh
40   USE dynspg_oce      ! pressure gradient schemes
41   USE dynspg_flt      ! pressure gradient schemes
42   USE dynspg_exp      ! pressure gradient schemes
43   USE dynspg_ts       ! pressure gradient schemes
44   USE sol_oce         ! ocean solver variables
45   USE lib_mpp         ! MPP library
46   USE wrk_nemo        ! Memory allocation
47   USE timing          ! Timing
48   USE asminc
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.3 , NEMO Consortium (2010)
60   !! $Id$
61   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
62   !!----------------------------------------------------------------------
63CONTAINS
64
65   SUBROUTINE istate_init
66      !!----------------------------------------------------------------------
67      !!                   ***  ROUTINE istate_init  ***
68      !!
69      !! ** Purpose :   Initialization of the dynamics and tracer fields.
70      !!----------------------------------------------------------------------
71      ! - ML - needed for initialization of e3t_b
72      INTEGER  ::  jk     ! dummy loop indice
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
84      rhd  (:,:,:  ) = 0.e0
85      rhop (:,:,:  ) = 0.e0
86      rn2  (:,:,:  ) = 0.e0
87      tsa  (:,:,:,:) = 0.e0
88      rn2b (:,:,:  ) = 0.e0
89
90      IF( ln_rstart ) THEN                    ! Restart from a file
91         !                                    ! -------------------
92         neuler = 1                              ! Set time-step indicator at nit000 (leap-frog)
93         CALL rst_read                           ! Read the restart file
94         !                                       ! define e3u_b, e3v_b from e3t_b read in restart file
95         CALL dom_vvl_2( nit000, fse3u_b(:,:,:), fse3v_b(:,:,:) )
96         CALL day_init                           ! model calendar (using both namelist and restart infos)
97      ELSE
98         !                                    ! Start from rest
99         !                                    ! ---------------
100         numror = 0                              ! define numror = 0 -> no restart file to read
101         neuler = 0                              ! Set time-step indicator at nit000 (euler forward)
102         CALL day_init                           ! model calendar (using both namelist and restart infos)
103         !                                       ! Initialization of ocean to zero
104         !   before fields      !       now fields
105         sshb (:,:)   = 0._wp   ;   sshn (:,:)   = 0._wp
106         ub   (:,:,:) = 0._wp   ;   un   (:,:,:) = 0._wp
107         vb   (:,:,:) = 0._wp   ;   vn   (:,:,:) = 0._wp
108         rotb (:,:,:) = 0._wp   ;   rotn (:,:,:) = 0._wp
109         hdivb(:,:,:) = 0._wp   ;   hdivn(:,:,:) = 0._wp
110         !
111         IF( cp_cfg == 'eel' ) THEN
112            CALL istate_eel                      ! EEL   configuration : start from pre-defined U,V T-S fields
113         ELSEIF( cp_cfg == 'gyre' ) THEN
114            CALL istate_gyre                     ! GYRE  configuration : start from pre-defined T-S fields
115         ELSEIF( cp_cfg == 'seabass' ) THEN
116            CALL istate_seabass
117         ELSEIF( ln_tsd_init      ) THEN         ! Initial T-S fields read in files
118            CALL dta_tsd( nit000, tsb )                  ! read 3D T and S data at nit000
119            tsn(:,:,:,:) = tsb(:,:,:,:)
120            !
121         ELSE                                    ! Initial T-S fields defined analytically
122            CALL istate_t_s
123         ENDIF
124         !
125         CALL eos( tsb, rhd, rhop )        ! before potential and in situ densities
126#if ! defined key_c1d
127         IF( ln_zps )   CALL zps_hde( nit000, jpts, tsb, gtsu, gtsv,  & ! zps: before hor. gradient
128            &                                       rhd, gru , grv  )   ! of t,s,rd at ocean bottom
129#endif
130         !
131         ! - ML - sshn could be modified by istate_eel, so that initialization of fse3t_b is done here
132         IF( lk_vvl ) THEN
133            DO jk = 1, jpk
134               fse3t_b(:,:,jk) = fse3t_n(:,:,jk)
135            ENDDO
136         ENDIF
137         !                                       ! define e3u_b, e3v_b from e3t_b initialized in domzgr
138         CALL dom_vvl_2( nit000, fse3u_b(:,:,:), fse3v_b(:,:,:) )
139         !
140      ENDIF
141      !
142      IF( lk_agrif ) THEN                  ! read free surface arrays in restart file
143         IF( ln_rstart ) THEN
144            IF( lk_dynspg_flt )  THEN      ! read or initialize the following fields
145               !                           ! gcx, gcxb for agrif_opa_init
146               IF( sol_oce_alloc()  > 0 )   CALL ctl_stop('agrif sol_oce_alloc: allocation of arrays failed')
147               CALL flt_rst( nit000, 'READ' )
148            ENDIF
149         ENDIF                             ! explicit case not coded yet with AGRIF
150      ENDIF
151      !
152      IF ( lk_asminc .AND. ln_asmdin ) THEN
153         neuler = 0            ! restart with an euler from the corrected background
154      ENDIF
155
156      IF( nn_timing == 1 )  CALL timing_stop('istate_init')
157      !
158   END SUBROUTINE istate_init
159
160   SUBROUTINE istate_t_s
161      !!---------------------------------------------------------------------
162      !!                  ***  ROUTINE istate_t_s  ***
163      !!
164      !! ** Purpose :   Intialization of the temperature field with an
165      !!      analytical profile or a file (i.e. in EEL configuration)
166      !!
167      !! ** Method  : - temperature: use Philander analytic profile
168      !!              - salinity   : use to a constant value 35.5
169      !!
170      !! References :  Philander ???
171      !!----------------------------------------------------------------------
172      INTEGER  :: ji, jj, jk
173      REAL(wp) ::   zsal = 35.50
174      !!----------------------------------------------------------------------
175      !
176      IF(lwp) WRITE(numout,*)
177      IF(lwp) WRITE(numout,*) 'istate_t_s : Philander s initial temperature profile'
178      IF(lwp) WRITE(numout,*) '~~~~~~~~~~   and constant salinity (',zsal,' psu)'
179      !
180      DO jk = 1, jpk
181         tsn(:,:,jk,jp_tem) = (  ( ( 7.5 - 0. * ABS( gphit(:,:) )/30. ) * ( 1.-TANH((fsdept(:,:,jk)-80.)/30.) )   &
182            &                + 10. * ( 5000. - fsdept(:,:,jk) ) /5000.)  ) * tmask(:,:,jk)
183         tsb(:,:,jk,jp_tem) = tsn(:,:,jk,jp_tem)
184      END DO
185      tsn(:,:,:,jp_sal) = zsal * tmask(:,:,:)
186      tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)
187      !
188   END SUBROUTINE istate_t_s
189
190
191   SUBROUTINE istate_eel
192      !!----------------------------------------------------------------------
193      !!                   ***  ROUTINE istate_eel  ***
194      !!
195      !! ** Purpose :   Initialization of the dynamics and tracers for EEL R5
196      !!      configuration (channel with or without a topographic bump)
197      !!
198      !! ** Method  : - set temprature field
199      !!              - set salinity field
200      !!              - set velocity field including horizontal divergence
201      !!                and relative vorticity fields
202      !!----------------------------------------------------------------------
203      USE divcur     ! hor. divergence & rel. vorticity      (div_cur routine)
204      USE iom
205
206      INTEGER  ::   inum              ! temporary logical unit
207      INTEGER  ::   ji, jj, jk        ! dummy loop indices
208      INTEGER  ::   ijloc
209      REAL(wp) ::   zh1, zh2, zslope, zcst, zfcor   ! temporary scalars
210      REAL(wp) ::   zt1  = 15._wp                   ! surface temperature value (EEL R5)
211      REAL(wp) ::   zt2  =  5._wp                   ! bottom  temperature value (EEL R5)
212      REAL(wp) ::   zsal = 35.0_wp                  ! constant salinity (EEL R2, R5 and R6)
213      REAL(wp) ::   zueel = 0.1_wp                  ! constant uniform zonal velocity (EEL R5)
214      REAL(wp), DIMENSION(jpiglo,jpjglo) ::   zssh  ! initial ssh over the global domain
215      !!----------------------------------------------------------------------
216
217      SELECT CASE ( jp_cfg )
218         !                                              ! ====================
219         CASE ( 5 )                                     ! EEL R5 configuration
220            !                                           ! ====================
221            !
222            ! set temperature field with a linear profile
223            ! -------------------------------------------
224            IF(lwp) WRITE(numout,*)
225            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: linear temperature profile'
226            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
227            !
228            zh1 = gdept_0(  1  )
229            zh2 = gdept_0(jpkm1)
230            !
231            zslope = ( zt1 - zt2 ) / ( zh1 - zh2 )
232            zcst   = ( zt1 * ( zh1 - zh2) - ( zt1 - zt2 ) * zh1 ) / ( zh1 - zh2 )
233            !
234            DO jk = 1, jpk
235               tsn(:,:,jk,jp_tem) = ( zt2 + zt1 * exp( - fsdept(:,:,jk) / 1000 ) ) * tmask(:,:,jk)
236               tsb(:,:,jk,jp_tem) = tsn(:,:,jk,jp_tem)
237            END DO
238            !
239            IF(lwp) CALL prizre( tsn(:,:,:,jp_tem), jpi   , jpj   , jpk   , jpj/2 ,   &
240               &                             1     , jpi   , 5     , 1     , jpk   ,   &
241               &                             1     , 1.    , numout                  )
242            !
243            ! set salinity field to a constant value
244            ! --------------------------------------
245            IF(lwp) WRITE(numout,*)
246            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal
247            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
248            !
249            tsn(:,:,:,jp_sal) = zsal * tmask(:,:,:)
250            tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)
251            !
252            ! set the dynamics: U,V, hdiv, rot (and ssh if necessary)
253            ! ----------------
254            ! Start EEL5 configuration with barotropic geostrophic velocities
255            ! according the sshb and sshn SSH imposed.
256            ! we assume a uniform grid (hence the use of e1t(1,1) for delta_y)
257            ! we use the Coriolis frequency at mid-channel.
258            ub(:,:,:) = zueel * umask(:,:,:)
259            un(:,:,:) = ub(:,:,:)
260            ijloc = mj0(INT(jpjglo-1)/2)
261            zfcor = ff(1,ijloc)
262            !
263            DO jj = 1, jpjglo
264               zssh(:,jj) = - (FLOAT(jj)- FLOAT(jpjglo-1)/2.)*zueel*e1t(1,1)*zfcor/grav
265            END DO
266            !
267            IF(lwp) THEN
268               WRITE(numout,*) ' Uniform zonal velocity for EEL R5:',zueel
269               WRITE(numout,*) ' Geostrophic SSH profile as a function of y:'
270               WRITE(numout,'(12(1x,f6.2))') zssh(1,:)
271            ENDIF
272            !
273            DO jj = 1, nlcj
274               DO ji = 1, nlci
275                  sshb(ji,jj) = zssh( mig(ji) , mjg(jj) ) * tmask(ji,jj,1)
276               END DO
277            END DO
278            sshb(nlci+1:jpi,      :   ) = 0.e0      ! set to zero extra mpp columns
279            sshb(      :   ,nlcj+1:jpj) = 0.e0      ! set to zero extra mpp rows
280            !
281            sshn(:,:) = sshb(:,:)                   ! set now ssh to the before value
282            !
283            IF( nn_rstssh /= 0 ) THEN
284               nn_rstssh = 0                        ! hand-made initilization of ssh
285               CALL ctl_warn( 'istate_eel: force nn_rstssh = 0' )
286            ENDIF
287            !
288            CALL div_cur( nit000 )                  ! horizontal divergence and relative vorticity (curl)
289            ! N.B. the vertical velocity will be computed from the horizontal divergence field
290            ! in istate by a call to wzv routine
291
292
293            !                                     ! ==========================
294         CASE ( 2 , 6 )                           ! EEL R2 or R6 configuration
295            !                                     ! ==========================
296            !
297            ! set temperature field with a NetCDF file
298            ! ----------------------------------------
299            IF(lwp) WRITE(numout,*)
300            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R2 or R6: read initial temperature in a NetCDF file'
301            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
302            !
303            CALL iom_open ( 'eel.initemp', inum )
304            CALL iom_get ( inum, jpdom_data, 'initemp', tsb(:,:,:,jp_tem) ) ! read before temprature (tb)
305            CALL iom_close( inum )
306            !
307            tsn(:,:,:,jp_tem) = tsb(:,:,:,jp_tem)                            ! set nox temperature to tb
308            !
309            IF(lwp) CALL prizre( tsn(:,:,:,jp_tem), jpi   , jpj   , jpk   , jpj/2 ,   &
310               &                            1     , jpi   , 5     , 1     , jpk   ,   &
311               &                            1     , 1.    , numout                  )
312            !
313            ! set salinity field to a constant value
314            ! --------------------------------------
315            IF(lwp) WRITE(numout,*)
316            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal
317            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
318            !
319            tsn(:,:,:,jp_sal) = zsal * tmask(:,:,:)
320            tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)
321            !
322            !                                    ! ===========================
323         CASE DEFAULT                            ! NONE existing configuration
324            !                                    ! ===========================
325            WRITE(ctmp1,*) 'EEL with a ', jp_cfg,' km resolution is not coded'
326            CALL ctl_stop( ctmp1 )
327            !
328      END SELECT
329      !
330   END SUBROUTINE istate_eel
331
332
333   SUBROUTINE istate_gyre
334      !!----------------------------------------------------------------------
335      !!                   ***  ROUTINE istate_gyre  ***
336      !!
337      !! ** Purpose :   Initialization of the dynamics and tracers for GYRE
338      !!      configuration (double gyre with rotated domain)
339      !!
340      !! ** Method  : - set temprature field
341      !!              - set salinity field
342      !!----------------------------------------------------------------------
343      INTEGER :: ji, jj, jk  ! dummy loop indices
344      INTEGER            ::   inum          ! temporary logical unit
345      INTEGER, PARAMETER ::   ntsinit = 0   ! (0/1) (analytical/input data files) T&S initialization
346      !!----------------------------------------------------------------------
347
348      SELECT CASE ( ntsinit)
349
350      CASE ( 0 )                  ! analytical T/S profil deduced from LEVITUS
351         IF(lwp) WRITE(numout,*)
352         IF(lwp) WRITE(numout,*) 'istate_gyre : initial analytical T and S profil deduced from LEVITUS '
353         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
354
355         DO jk = 1, jpk
356            DO jj = 1, jpj
357               DO ji = 1, jpi
358                  tsn(ji,jj,jk,jp_tem) = (  16. - 12. * TANH( (fsdept(ji,jj,jk) - 400) / 700 )         )   &
359                       &           * (-TANH( (500-fsdept(ji,jj,jk)) / 150 ) + 1) / 2               &
360                       &       + (      15. * ( 1. - TANH( (fsdept(ji,jj,jk)-50.) / 1500.) )       &
361                       &                - 1.4 * TANH((fsdept(ji,jj,jk)-100.) / 100.)               &
362                       &                + 7.  * (1500. - fsdept(ji,jj,jk)) / 1500.             )   &
363                       &           * (-TANH( (fsdept(ji,jj,jk) - 500) / 150) + 1) / 2
364                  tsn(ji,jj,jk,jp_tem) = tsn(ji,jj,jk,jp_tem) * tmask(ji,jj,jk)
365                  tsb(ji,jj,jk,jp_tem) = tsn(ji,jj,jk,jp_tem)
366
367                  tsn(ji,jj,jk,jp_sal) =  (  36.25 - 1.13 * TANH( (fsdept(ji,jj,jk) - 305) / 460 )  )  &
368                     &              * (-TANH((500 - fsdept(ji,jj,jk)) / 150) + 1) / 2          &
369                     &          + (  35.55 + 1.25 * (5000. - fsdept(ji,jj,jk)) / 5000.         &
370                     &                - 1.62 * TANH( (fsdept(ji,jj,jk) - 60.  ) / 650. )       &
371                     &                + 0.2  * TANH( (fsdept(ji,jj,jk) - 35.  ) / 100. )       &
372                     &                + 0.2  * TANH( (fsdept(ji,jj,jk) - 1000.) / 5000.)    )  &
373                     &              * (-TANH((fsdept(ji,jj,jk) - 500) / 150) + 1) / 2
374                  tsn(ji,jj,jk,jp_sal) = tsn(ji,jj,jk,jp_sal) * tmask(ji,jj,jk)
375                  tsb(ji,jj,jk,jp_sal) = tsn(ji,jj,jk,jp_sal)
376               END DO
377            END DO
378         END DO
379
380      CASE ( 1 )                  ! T/S data fields read in dta_tem.nc/data_sal.nc files
381         IF(lwp) WRITE(numout,*)
382         IF(lwp) WRITE(numout,*) 'istate_gyre : initial T and S read from dta_tem.nc/data_sal.nc files'
383         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
384         IF(lwp) WRITE(numout,*) '              NetCDF FORMAT'
385
386         ! Read temperature field
387         ! ----------------------
388         CALL iom_open ( 'data_tem', inum )
389         CALL iom_get ( inum, jpdom_data, 'votemper', tsn(:,:,:,jp_tem) )
390         CALL iom_close( inum )
391
392         tsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem) * tmask(:,:,:)
393         tsb(:,:,:,jp_tem) = tsn(:,:,:,jp_tem)
394
395         ! Read salinity field
396         ! -------------------
397         CALL iom_open ( 'data_sal', inum )
398         CALL iom_get ( inum, jpdom_data, 'vosaline', tsn(:,:,:,jp_sal) )
399         CALL iom_close( inum )
400
401         tsn(:,:,:,jp_sal) = tsn(:,:,:,jp_sal) * tmask(:,:,:)
402         tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)
403
404      END SELECT
405
406      IF(lwp) THEN
407         WRITE(numout,*)
408         WRITE(numout,*) '              Initial temperature and salinity profiles:'
409         WRITE(numout, "(9x,' level   gdept_0   temperature   salinity   ')" )
410         WRITE(numout, "(10x, i4, 3f10.2)" ) ( jk, gdept_0(jk), tsn(2,2,jk,jp_tem), tsn(2,2,jk,jp_sal), jk = 1, jpk )
411      ENDIF
412
413   END SUBROUTINE istate_gyre
414
415
416
417   SUBROUTINE istate_seabass
418      !!----------------------------------------------------------------------
419      !!                   ***  ROUTINE istate_seabass  ***
420      !!
421      !! ** Purpose :   Initialization of the dynamics and tracers for seabass
422      !!      configuration (double gyre)
423      !!
424      !! ** Method  : - set temperature field following Chassignet and Gent, JPO
425      !!                    21, pp1290-1299, 1991, and the law
426      !!                    rho/rho0=1-2.e-4(T-T0)
427      !!              - set salinity field constant
428      !!
429      !!----------------------------------------------------------------------
430      !! * Local variables
431      INTEGER :: ji, jj, jk     ! dummy loop indices
432      REAL(wp) ::   zsal = 35.5
433      !!----------------------------------------------------------------------
434
435      IF(lwp) WRITE(numout,*)
436      IF(lwp) WRITE(numout,*) 'istate_seabass : initial analytical T and constant S profiles '
437      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
438
439      DO jk = 1, jpk
440         DO jj = 1, jpj
441            DO ji = 1, jpi
442               tsn(ji,jj,jk,jp_tem) = ( 25.+5.9e-5*800./9.81/2.e-4*   &
443               &    (exp(-gdept_0(jk)/800.)-1.))  * tmask(ji,jj,jk)
444               tsb(ji,jj,jk,jp_tem) = tsn(ji,jj,jk, jp_tem)
445          END DO
446        END DO
447      END DO
448
449      tsn(:,:,:,jp_sal) = zsal  * tmask(:,:,:)
450      tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal)
451
452      IF(lwp) THEN
453         WRITE(numout,*)
454         WRITE(numout,*) '              Initial temperature and salinity profiles:'
455         WRITE(numout, "(9x,' level   gdept   temperature   salinity   ')" )
456         WRITE(numout, "(10x, i4, 3f10.2)" ) ( jk, gdept_0(jk), tsn(2,2,jk,jp_tem), tsn(2,2,jk,jp_sal), jk = 1, jpk )
457      ENDIF
458
459
460   END SUBROUTINE istate_seabass
461
462
463
464   SUBROUTINE istate_uvg
465      !!----------------------------------------------------------------------
466      !!                  ***  ROUTINE istate_uvg  ***
467      !!
468      !! ** Purpose :   Compute the geostrophic velocities from (tn,sn) fields
469      !!
470      !! ** Method  :   Using the hydrostatic hypothesis the now hydrostatic
471      !!      pressure is computed by integrating the in-situ density from the
472      !!      surface to the bottom.
473      !!                 p=integral [ rau*g dz ]
474      !!----------------------------------------------------------------------
475      USE dynspg          ! surface pressure gradient             (dyn_spg routine)
476      USE divcur          ! hor. divergence & rel. vorticity      (div_cur routine)
477      USE lbclnk          ! ocean lateral boundary condition (or mpp link)
478
479      INTEGER ::   ji, jj, jk        ! dummy loop indices
480      INTEGER ::   indic             ! ???
481      REAL(wp) ::   zmsv, zphv, zmsu, zphu, zalfg     ! temporary scalars
482      REAL(wp), POINTER, DIMENSION(:,:,:) :: zprn
483      !!----------------------------------------------------------------------
484      !
485      CALL wrk_alloc( jpi, jpj, jpk, zprn)
486      !
487      IF(lwp) WRITE(numout,*)
488      IF(lwp) WRITE(numout,*) 'istate_uvg : Start from Geostrophy'
489      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
490
491      ! Compute the now hydrostatic pressure
492      ! ------------------------------------
493
494      zalfg = 0.5 * grav * rau0
495
496      zprn(:,:,1) = zalfg * fse3w(:,:,1) * ( 1 + rhd(:,:,1) )       ! Surface value
497
498      DO jk = 2, jpkm1                                              ! Vertical integration from the surface
499         zprn(:,:,jk) = zprn(:,:,jk-1)   &
500            &         + zalfg * fse3w(:,:,jk) * ( 2. + rhd(:,:,jk) + rhd(:,:,jk-1) )
501      END DO
502
503      ! Compute geostrophic balance
504      ! ---------------------------
505      DO jk = 1, jpkm1
506         DO jj = 2, jpjm1
507            DO ji = fs_2, fs_jpim1   ! vertor opt.
508               zmsv = 1. / MAX(  umask(ji-1,jj+1,jk) + umask(ji  ,jj+1,jk)   &
509                               + umask(ji-1,jj  ,jk) + umask(ji  ,jj  ,jk) , 1.  )
510               zphv = ( zprn(ji  ,jj+1,jk) - zprn(ji-1,jj+1,jk) ) * umask(ji-1,jj+1,jk) / e1u(ji-1,jj+1)   &
511                    + ( zprn(ji+1,jj+1,jk) - zprn(ji  ,jj+1,jk) ) * umask(ji  ,jj+1,jk) / e1u(ji  ,jj+1)   &
512                    + ( zprn(ji  ,jj  ,jk) - zprn(ji-1,jj  ,jk) ) * umask(ji-1,jj  ,jk) / e1u(ji-1,jj  )   &
513                    + ( zprn(ji+1,jj  ,jk) - zprn(ji  ,jj  ,jk) ) * umask(ji  ,jj  ,jk) / e1u(ji  ,jj  )
514               zphv = 1. / rau0 * zphv * zmsv * vmask(ji,jj,jk)
515
516               zmsu = 1. / MAX(  vmask(ji+1,jj  ,jk) + vmask(ji  ,jj  ,jk)   &
517                               + vmask(ji+1,jj-1,jk) + vmask(ji  ,jj-1,jk) , 1.  )
518               zphu = ( zprn(ji+1,jj+1,jk) - zprn(ji+1,jj  ,jk) ) * vmask(ji+1,jj  ,jk) / e2v(ji+1,jj  )   &
519                    + ( zprn(ji  ,jj+1,jk) - zprn(ji  ,jj  ,jk) ) * vmask(ji  ,jj  ,jk) / e2v(ji  ,jj  )   &
520                    + ( zprn(ji+1,jj  ,jk) - zprn(ji+1,jj-1,jk) ) * vmask(ji+1,jj-1,jk) / e2v(ji+1,jj-1)   &
521                    + ( zprn(ji  ,jj  ,jk) - zprn(ji  ,jj-1,jk) ) * vmask(ji  ,jj-1,jk) / e2v(ji  ,jj-1)
522               zphu = 1. / rau0 * zphu * zmsu * umask(ji,jj,jk)
523
524               ! Compute the geostrophic velocities
525               un(ji,jj,jk) = -2. * zphu / ( ff(ji,jj) + ff(ji  ,jj-1) )
526               vn(ji,jj,jk) =  2. * zphv / ( ff(ji,jj) + ff(ji-1,jj  ) )
527            END DO
528         END DO
529      END DO
530
531      IF(lwp) WRITE(numout,*) '         we force to zero bottom velocity'
532
533      ! Susbtract the bottom velocity (level jpk-1 for flat bottom case)
534      ! to have a zero bottom velocity
535
536      DO jk = 1, jpkm1
537         un(:,:,jk) = ( un(:,:,jk) - un(:,:,jpkm1) ) * umask(:,:,jk)
538         vn(:,:,jk) = ( vn(:,:,jk) - vn(:,:,jpkm1) ) * vmask(:,:,jk)
539      END DO
540
541      CALL lbc_lnk( un, 'U', -1. )
542      CALL lbc_lnk( vn, 'V', -1. )
543
544      ub(:,:,:) = un(:,:,:)
545      vb(:,:,:) = vn(:,:,:)
546
547      ! WARNING !!!!!
548      ! after initializing u and v, we need to calculate the initial streamfunction bsf.
549      ! Otherwise, only the trend will be computed and the model will blow up (inconsistency).
550      ! to do that, we call dyn_spg with a special trick:
551      ! we fill ua and va with the velocities divided by dt, and the streamfunction will be brought to the
552      ! right value assuming the velocities have been set up in one time step.
553      ! we then set bsfd to zero (first guess for next step is d(psi)/dt = 0.)
554      !  sets up s false trend to calculate the barotropic streamfunction.
555
556      ua(:,:,:) = ub(:,:,:) / rdt
557      va(:,:,:) = vb(:,:,:) / rdt
558
559      ! calls dyn_spg. we assume euler time step, starting from rest.
560      indic = 0
561      CALL dyn_spg( nit000, indic )       ! surface pressure gradient
562
563      ! the new velocity is ua*rdt
564
565      CALL lbc_lnk( ua, 'U', -1. )
566      CALL lbc_lnk( va, 'V', -1. )
567
568      ub(:,:,:) = ua(:,:,:) * rdt
569      vb(:,:,:) = va(:,:,:) * rdt
570      ua(:,:,:) = 0.e0
571      va(:,:,:) = 0.e0
572      un(:,:,:) = ub(:,:,:)
573      vn(:,:,:) = vb(:,:,:)
574
575      ! Compute the divergence and curl
576
577      CALL div_cur( nit000 )            ! now horizontal divergence and curl
578
579      hdivb(:,:,:) = hdivn(:,:,:)       ! set the before to the now value
580      rotb (:,:,:) = rotn (:,:,:)       ! set the before to the now value
581      !
582      CALL wrk_dealloc( jpi, jpj, jpk, zprn)
583      !
584   END SUBROUTINE istate_uvg
585
586   !!=====================================================================
587END MODULE istate
Note: See TracBrowser for help on using the repository browser.