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

source: branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/DOM/istate.F90 @ 4400

Last change on this file since 4400 was 3837, checked in by trackstand2, 11 years ago

Merge of finiss

  • Property svn:keywords set to Id
File size: 28.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   !!----------------------------------------------------------------------
16
17   !!----------------------------------------------------------------------
18   !!   istate_init   : initial state setting
19   !!   istate_tem    : analytical profile for initial Temperature
20   !!   istate_sal    : analytical profile for initial Salinity
21   !!   istate_eel    : initial state setting of EEL R5 configuration
22   !!   istate_gyre   : initial state setting of GYRE configuration
23   !!   istate_uvg    : initial velocity in geostropic balance
24   !!----------------------------------------------------------------------
25   USE oce             ! ocean dynamics and active tracers
26   USE dom_oce         ! ocean space and time domain
27   USE daymod          ! calendar
28   USE eosbn2          ! eq. of state, Brunt Vaisala frequency (eos     routine)
29   USE ldftra_oce      ! ocean active tracers: lateral physics
30   USE zdf_oce         ! ocean vertical physics
31   USE phycst          ! physical constants
32   USE dtatem          ! temperature data                 (dta_tem routine)
33   USE dtasal          ! salinity data                    (dta_sal 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 traswp          ! Swap arrays                      (tra_swp routine)
45   USE lib_mpp         ! MPP library
46
47   IMPLICIT NONE
48   PRIVATE
49
50   PUBLIC   istate_init   ! routine called by step.F90
51
52   !! * Control permutation of array indices
53#  include "oce_ftrans.h90"
54#  include "dom_oce_ftrans.h90"
55#  include "ldftra_oce_ftrans.h90"
56#  include "zdf_oce_ftrans.h90"
57#  include "dtatem_ftrans.h90"
58#  include "dtasal_ftrans.h90"
59#  include "domvvl_ftrans.h90"
60
61   !! * Substitutions
62#  include "domzgr_substitute.h90"
63#  include "vectopt_loop_substitute.h90"
64   !!----------------------------------------------------------------------
65   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
66   !! $Id$
67   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
68   !!----------------------------------------------------------------------
69CONTAINS
70
71   SUBROUTINE istate_init
72      !!----------------------------------------------------------------------
73      !!                   ***  ROUTINE istate_init  ***
74      !!
75      !! ** Purpose :   Initialization of the dynamics and tracer fields.
76      !!----------------------------------------------------------------------
77      ! - ML - needed for initialization of e3t_b
78      INTEGER  ::  ji, jj, jk     ! dummy loop indices
79
80      IF(lwp) WRITE(numout,*)
81      IF(lwp) WRITE(numout,*) 'istate_ini : Initialization of the dynamics and tracers'
82      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
83
84      rhd  (:,:,:) = 0.e0
85      rhop (:,:,:) = 0.e0
86      rn2  (:,:,:) = 0.e0 
87      ta   (:,:,:) = 0.e0   
88      sa   (:,:,:) = 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         CALL tra_swap                           ! swap 3D arrays (t,s)  in a 4D array (ts)
95         CALL day_init                           ! model calendar (using both namelist and restart infos)
96      ELSE
97         !                                    ! Start from rest
98         !                                    ! ---------------
99         numror = 0                              ! define numror = 0 -> no restart file to read
100         neuler = 0                              ! Set time-step indicator at nit000 (euler forward)
101         CALL day_init                           ! model calendar (using both namelist and restart infos)
102         !                                       ! Initialization of ocean to zero
103         !   before fields     !       now fields     
104         sshb (:,:)   = 0.e0   ;   sshn (:,:)   = 0.e0
105         ub   (:,:,:) = 0.e0   ;   un   (:,:,:) = 0.e0
106         vb   (:,:,:) = 0.e0   ;   vn   (:,:,:) = 0.e0 
107         rotb (:,:,:) = 0.e0   ;   rotn (:,:,:) = 0.e0
108         hdivb(:,:,:) = 0.e0   ;   hdivn(:,:,:) = 0.e0
109         !
110         IF( cp_cfg == 'eel' ) THEN
111            CALL istate_eel                      ! EEL   configuration : start from pre-defined U,V T-S fields
112         ELSEIF( cp_cfg == 'gyre' ) THEN         
113            CALL istate_gyre                     ! GYRE  configuration : start from pre-defined T-S fields
114         ELSE
115            !                                    ! Other configurations: Initial T-S fields
116#if defined key_dtatem
117            CALL dta_tem( nit000 )                  ! read 3D temperature data
118            tb(:,:,:) = t_dta(:,:,:)   ;   tn(:,:,:) = t_dta(:,:,:)
119           
120#else
121            IF(lwp) WRITE(numout,*)                 ! analytical temperature profile
122            IF(lwp) WRITE(numout,*)'             Temperature initialization using an analytic profile'
123            CALL istate_tem
124#endif
125#if defined key_dtasal
126            CALL dta_sal( nit000 )                  ! read 3D salinity data
127            sb(:,:,:) = s_dta(:,:,:)   ;   sn(:,:,:) = s_dta(:,:,:)
128#else
129            ! No salinity data
130            IF(lwp)WRITE(numout,*)                  ! analytical salinity profile
131            IF(lwp)WRITE(numout,*)'             Salinity initialisation using a constant value'
132            CALL istate_sal
133#endif
134         ENDIF
135         !
136         CALL tra_swap                     ! swap 3D arrays (tb,sb,tn,sn)  in a 4D array
137         CALL eos( tsb, rhd, rhop )        ! before potential and in situ densities
138#if ! defined key_c1d
139         IF( ln_zps )   CALL zps_hde( nit000, jpts, tsb, gtsu, gtsv,  & ! zps: before hor. gradient
140            &                                       rhd, gru , grv  )   ! of t,s,rd at ocean bottom
141#endif
142         !   
143         ! - ML - sshn could be modified by istate_eel, so that initialization of fse3t_b is done here
144         IF( lk_vvl ) THEN
145#if defined key_z_first
146            ! DCSE_NEMO: can't use implicit loop over k here because the domzgr_substitute.h90
147            ! file causes the line below to be expanded to:
148            !     e3t_b(:,:,:) = (e3t(:,:,:)*(1.+sshn(:,:)*mut(:,:,:)))
149            ! which contains non-conforming array expressions.
150            DO jk = 1, jpk
151               fse3t_b(:,:,jk) = fse3t_n(:,:,jk)
152            END DO
153#else
154            DO jk = 1, jpk
155               fse3t_b(:,:,jk) = fse3t_n(:,:,jk)
156            ENDDO
157#endif
158         ENDIF
159         !
160      ENDIF
161      !
162      IF( lk_agrif ) THEN                  ! read free surface arrays in restart file
163         IF( ln_rstart ) THEN
164            IF( lk_dynspg_flt )   CALL flt_rst( nit000, 'READ' )      ! read or initialize the following fields
165            !                                                         ! gcx, gcxb for agrif_opa_init
166         ENDIF                                                        ! explicit case not coded yet with AGRIF
167      ENDIF
168
169   END SUBROUTINE istate_init
170
171
172   SUBROUTINE istate_tem
173      !!---------------------------------------------------------------------
174      !!                  ***  ROUTINE istate_tem  ***
175      !!   
176      !! ** Purpose :   Intialization of the temperature field with an
177      !!      analytical profile or a file (i.e. in EEL configuration)
178      !!
179      !! ** Method  :   Use Philander analytic profile of temperature
180      !!
181      !! References :  Philander ???
182      !!----------------------------------------------------------------------
183      INTEGER :: ji, jj, jk
184      !!----------------------------------------------------------------------
185      !
186      IF(lwp) WRITE(numout,*)
187      IF(lwp) WRITE(numout,*) 'istate_tem : initial temperature profile'
188      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
189      !
190#if defined key_z_first
191      DO jj = 1, jpj
192         DO ji = 1, jpi
193            DO jk = 1, jpk
194#else
195      DO jk = 1, jpk
196         DO jj = 1, jpj
197            DO ji = 1, jpi
198#endif
199               tn(ji,jj,jk) = (  ( ( 7.5 - 0.*ABS(gphit(ji,jj))/30. )   &
200                  &               *( 1.-TANH((fsdept(ji,jj,jk)-80.)/30.) )   &
201                  &            + 10.*(5000.-fsdept(ji,jj,jk))/5000.)  ) * tmask(ji,jj,jk)
202               tb(ji,jj,jk) = tn(ji,jj,jk)
203          END DO
204        END DO
205      END DO
206      !
207      IF(lwp) CALL prizre( tn    , jpi   , jpj   , jpk   , jpj/2 ,   &
208         &                 1     , jpi   , 5     , 1     , jpk   ,   &
209         &                 1     , 1.    , numout                  )
210      !
211   END SUBROUTINE istate_tem
212
213
214   SUBROUTINE istate_sal
215      !!---------------------------------------------------------------------
216      !!                  ***  ROUTINE istate_sal  ***
217      !!
218      !! ** Purpose :   Intialize the salinity field with an analytic profile
219      !!
220      !! ** Method  :   Use to a constant value 35.5
221      !!             
222      !! ** Action  :   Initialize sn and sb
223      !!----------------------------------------------------------------------
224      REAL(wp) ::   zsal = 35.50_wp
225      !!----------------------------------------------------------------------
226      !
227      IF(lwp) WRITE(numout,*)
228      IF(lwp) WRITE(numout,*) 'istate_sal : initial salinity : ', zsal
229      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
230      !
231      sn(:,:,:) = zsal * tmask(:,:,:)
232      sb(:,:,:) = sn(:,:,:)
233      !
234   END SUBROUTINE istate_sal
235
236
237   SUBROUTINE istate_eel
238      !!----------------------------------------------------------------------
239      !!                   ***  ROUTINE istate_eel  ***
240      !!
241      !! ** Purpose :   Initialization of the dynamics and tracers for EEL R5
242      !!      configuration (channel with or without a topographic bump)
243      !!
244      !! ** Method  : - set temprature field
245      !!              - set salinity field
246      !!              - set velocity field including horizontal divergence
247      !!                and relative vorticity fields
248      !!----------------------------------------------------------------------
249      USE divcur     ! hor. divergence & rel. vorticity      (div_cur routine)
250      USE iom
251 
252      INTEGER  ::   inum              ! temporary logical unit
253      INTEGER  ::   ji, jj, jk        ! dummy loop indices
254      INTEGER  ::   ijloc
255      REAL(wp) ::   zh1, zh2, zslope, zcst, zfcor   ! temporary scalars
256      REAL(wp) ::   zt1  = 15._wp                   ! surface temperature value (EEL R5)
257      REAL(wp) ::   zt2  =  5._wp                   ! bottom  temperature value (EEL R5)
258      REAL(wp) ::   zsal = 35.0_wp                  ! constant salinity (EEL R2, R5 and R6)
259      REAL(wp) ::   zueel = 0.1_wp                  ! constant uniform zonal velocity (EEL R5)
260      REAL(wp), DIMENSION(jpiglo,jpjglo) ::   zssh  ! initial ssh over the global domain
261      !!----------------------------------------------------------------------
262
263      SELECT CASE ( jp_cfg ) 
264         !                                              ! ====================
265         CASE ( 5 )                                     ! EEL R5 configuration
266            !                                           ! ====================
267            !
268            ! set temperature field with a linear profile
269            ! -------------------------------------------
270            IF(lwp) WRITE(numout,*)
271            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: linear temperature profile'
272            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
273            !
274            zh1 = gdept_0(  1  )
275            zh2 = gdept_0(jpkm1)
276            !
277            zslope = ( zt1 - zt2 ) / ( zh1 - zh2 )
278            zcst   = ( zt1 * ( zh1 - zh2) - ( zt1 - zt2 ) * zh1 ) / ( zh1 - zh2 )
279            !
280#if defined key_z_first
281            DO jj = 1, jpj
282               DO ji = 1, jpi
283                  DO jk = 1, jpk
284                     tn(ji,jj,jk) = ( zt2 + zt1 * exp( - fsdept(ji,jj,jk) / 1000 ) ) * tmask(ji,jj,jk)
285                     tb(ji,jj,jk) = tn(ji,jj,jk)
286                  END DO
287               END DO
288            END DO
289#else
290            DO jk = 1, jpk
291               tn(:,:,jk) = ( zt2 + zt1 * exp( - fsdept(:,:,jk) / 1000 ) ) * tmask(:,:,jk)
292               tb(:,:,jk) = tn(:,:,jk)
293            END DO
294#endif
295            !
296            IF(lwp) CALL prizre( tn    , jpi   , jpj   , jpk   , jpj/2 ,   &
297               &                 1     , jpi   , 5     , 1     , jpk   ,   &
298               &                 1     , 1.    , numout                  )
299            !
300            ! set salinity field to a constant value
301            ! --------------------------------------
302            IF(lwp) WRITE(numout,*)
303            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal
304            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
305            !
306            sn(:,:,:) = zsal * tmask(:,:,:)
307            sb(:,:,:) = sn(:,:,:)
308            !
309            ! set the dynamics: U,V, hdiv, rot (and ssh if necessary)
310            ! ----------------
311            ! Start EEL5 configuration with barotropic geostrophic velocities
312            ! according the sshb and sshn SSH imposed.
313            ! we assume a uniform grid (hence the use of e1t(1,1) for delta_y)
314            ! we use the Coriolis frequency at mid-channel.   
315            ub(:,:,:) = zueel * umask(:,:,:)
316            un(:,:,:) = ub(:,:,:)
317            ijloc = mj0(INT(jpjglo-1)/2)
318            zfcor = ff(1,ijloc)
319            !
320            DO jj = 1, jpjglo
321               zssh(:,jj) = - (FLOAT(jj)- FLOAT(jpjglo-1)/2.)*zueel*e1t(1,1)*zfcor/grav 
322            END DO
323            !
324            IF(lwp) THEN
325               WRITE(numout,*) ' Uniform zonal velocity for EEL R5:',zueel
326               WRITE(numout,*) ' Geostrophic SSH profile as a function of y:'
327               WRITE(numout,'(12(1x,f6.2))') zssh(1,:)
328            ENDIF
329            !
330            DO jj = 1, nlcj
331               DO ji = 1, nlci
332#if defined key_z_first
333                  sshb(ji,jj) = zssh( mig(ji) , mjg(jj) ) * tmask_1(ji,jj)
334#else
335                  sshb(ji,jj) = zssh( mig(ji) , mjg(jj) ) * tmask(ji,jj,1)
336#endif
337               END DO
338            END DO
339            sshb(nlci+1:jpi,      :   ) = 0.e0      ! set to zero extra mpp columns
340            sshb(      :   ,nlcj+1:jpj) = 0.e0      ! set to zero extra mpp rows
341            !
342            sshn(:,:) = sshb(:,:)                   ! set now ssh to the before value
343            !
344            IF( nn_rstssh /= 0 ) THEN 
345               nn_rstssh = 0                        ! hand-made initilization of ssh
346               CALL ctl_warn( 'istate_eel: force nn_rstssh = 0' )
347            ENDIF
348            !
349            CALL div_cur( nit000 )                  ! horizontal divergence and relative vorticity (curl)
350            ! N.B. the vertical velocity will be computed from the horizontal divergence field
351            ! in istate by a call to wzv routine
352
353
354            !                                     ! ==========================
355         CASE ( 2 , 6 )                           ! EEL R2 or R6 configuration
356            !                                     ! ==========================
357            !
358            ! set temperature field with a NetCDF file
359            ! ----------------------------------------
360            IF(lwp) WRITE(numout,*)
361            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R2 or R6: read initial temperature in a NetCDF file'
362            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
363            !
364            CALL iom_open ( 'eel.initemp', inum )
365            CALL iom_get ( inum, jpdom_data, 'initemp', tb ) ! read before temprature (tb)
366            CALL iom_close( inum )
367            !
368            tn(:,:,:) = tb(:,:,:)                            ! set nox temperature to tb
369            !
370            IF(lwp) CALL prizre( tn    , jpi   , jpj   , jpk   , jpj/2 ,   &
371               &                 1     , jpi   , 5     , 1     , jpk   ,   &
372               &                 1     , 1.    , numout                  )
373            !
374            ! set salinity field to a constant value
375            ! --------------------------------------
376            IF(lwp) WRITE(numout,*)
377            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal
378            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
379            !
380            sn(:,:,:) = zsal * tmask(:,:,:)
381            sb(:,:,:) = sn(:,:,:)
382            !
383            !                                    ! ===========================
384         CASE DEFAULT                            ! NONE existing configuration
385            !                                    ! ===========================
386            WRITE(ctmp1,*) 'EEL with a ', jp_cfg,' km resolution is not coded'
387            CALL ctl_stop( ctmp1 )
388            !
389      END SELECT
390      !
391   END SUBROUTINE istate_eel
392
393
394   SUBROUTINE istate_gyre
395      !!----------------------------------------------------------------------
396      !!                   ***  ROUTINE istate_gyre  ***
397      !!
398      !! ** Purpose :   Initialization of the dynamics and tracers for GYRE
399      !!      configuration (double gyre with rotated domain)
400      !!
401      !! ** Method  : - set temprature field
402      !!              - set salinity field
403      !!----------------------------------------------------------------------
404      INTEGER :: ji, jj, jk  ! dummy loop indices
405      INTEGER            ::   inum          ! temporary logical unit
406      INTEGER, PARAMETER ::   ntsinit = 0   ! (0/1) (analytical/input data files) T&S initialization
407      !!----------------------------------------------------------------------
408
409      SELECT CASE ( ntsinit)
410
411      CASE ( 0 )                  ! analytical T/S profil deduced from LEVITUS
412         IF(lwp) WRITE(numout,*)
413         IF(lwp) WRITE(numout,*) 'istate_gyre : initial analytical T and S profil deduced from LEVITUS '
414         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
415
416#if defined key_z_first
417         DO jj = 1, jpj
418            DO ji = 1, jpi
419               DO jk = 1, jpk
420#else
421         DO jk = 1, jpk
422            DO jj = 1, jpj
423               DO ji = 1, jpi
424#endif
425                  tn(ji,jj,jk) = (  16. - 12. * TANH( (fsdept(ji,jj,jk) - 400) / 700 )         )   &
426                       &           * (-TANH( (500-fsdept(ji,jj,jk)) / 150 ) + 1) / 2               &
427                       &       + (      15. * ( 1. - TANH( (fsdept(ji,jj,jk)-50.) / 1500.) )       &
428                       &                - 1.4 * TANH((fsdept(ji,jj,jk)-100.) / 100.)               &   
429                       &                + 7.  * (1500. - fsdept(ji,jj,jk)) / 1500.             )   & 
430                       &           * (-TANH( (fsdept(ji,jj,jk) - 500) / 150) + 1) / 2
431                  tn(ji,jj,jk) = tn(ji,jj,jk) * tmask(ji,jj,jk)
432                  tb(ji,jj,jk) = tn(ji,jj,jk)
433
434                  sn(ji,jj,jk) =  (  36.25 - 1.13 * TANH( (fsdept(ji,jj,jk) - 305) / 460 )  )  &
435                     &              * (-TANH((500 - fsdept(ji,jj,jk)) / 150) + 1) / 2          &
436                     &          + (  35.55 + 1.25 * (5000. - fsdept(ji,jj,jk)) / 5000.         &
437                     &                - 1.62 * TANH( (fsdept(ji,jj,jk) - 60.  ) / 650. )       &
438                     &                + 0.2  * TANH( (fsdept(ji,jj,jk) - 35.  ) / 100. )       &
439                     &                + 0.2  * TANH( (fsdept(ji,jj,jk) - 1000.) / 5000.)    )  &
440                     &              * (-TANH((fsdept(ji,jj,jk) - 500) / 150) + 1) / 2 
441                  sn(ji,jj,jk) = sn(ji,jj,jk) * tmask(ji,jj,jk)
442                  sb(ji,jj,jk) = sn(ji,jj,jk)
443               END DO
444            END DO
445         END DO
446
447      CASE ( 1 )                  ! T/S data fields read in dta_tem.nc/data_sal.nc files
448         IF(lwp) WRITE(numout,*)
449         IF(lwp) WRITE(numout,*) 'istate_gyre : initial T and S read from dta_tem.nc/data_sal.nc files'
450         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
451         IF(lwp) WRITE(numout,*) '              NetCDF FORMAT'
452
453         ! Read temperature field
454         ! ----------------------
455         CALL iom_open ( 'data_tem', inum )
456         CALL iom_get ( inum, jpdom_data, 'votemper', tn ) 
457         CALL iom_close( inum )
458
459         tn(:,:,:) = tn(:,:,:) * tmask(:,:,:) 
460         tb(:,:,:) = tn(:,:,:)
461
462         ! Read salinity field
463         ! -------------------
464         CALL iom_open ( 'data_sal', inum )
465         CALL iom_get ( inum, jpdom_data, 'vosaline', sn ) 
466         CALL iom_close( inum )
467
468         sn(:,:,:)  = sn(:,:,:) * tmask(:,:,:) 
469         sb(:,:,:)  = sn(:,:,:)
470
471      END SELECT
472
473      IF(lwp) THEN
474         WRITE(numout,*)
475         WRITE(numout,*) '              Initial temperature and salinity profiles:'
476         WRITE(numout, "(9x,' level   gdept_0   temperature   salinity   ')" )
477         WRITE(numout, "(10x, i4, 3f10.2)" ) ( jk, gdept_0(jk), tn(2,2,jk), sn(2,2,jk), jk = 1, jpk )
478      ENDIF
479
480   END SUBROUTINE istate_gyre
481
482
483   SUBROUTINE istate_uvg
484      !!----------------------------------------------------------------------
485      !!                  ***  ROUTINE istate_uvg  ***
486      !!
487      !! ** Purpose :   Compute the geostrophic velocities from (tn,sn) fields
488      !!
489      !! ** Method  :   Using the hydrostatic hypothesis the now hydrostatic
490      !!      pressure is computed by integrating the in-situ density from the
491      !!      surface to the bottom.
492      !!                 p=integral [ rau*g dz ]
493      !!----------------------------------------------------------------------
494      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
495      USE wrk_nemo, ONLY:   zprn => wrk_3d_1    ! 3D workspace
496      !! DCSE_NEMO: wrk_3d_1 renamed, need additional directive
497!FTRANS zprn :I :I :z
498
499      USE dynspg          ! surface pressure gradient             (dyn_spg routine)
500      USE divcur          ! hor. divergence & rel. vorticity      (div_cur routine)
501      USE lbclnk          ! ocean lateral boundary condition (or mpp link)
502
503      INTEGER ::   ji, jj, jk        ! dummy loop indices
504      INTEGER ::   indic             ! ???
505      REAL(wp) ::   zmsv, zphv, zmsu, zphu, zalfg     ! temporary scalars
506      !!----------------------------------------------------------------------
507
508      IF(wrk_in_use(3, 1) ) THEN
509         CALL ctl_stop('istate_uvg: requested workspace array unavailable')   ;   RETURN
510      ENDIF
511
512      IF(lwp) WRITE(numout,*) 
513      IF(lwp) WRITE(numout,*) 'istate_uvg : Start from Geostrophy'
514      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
515
516      ! Compute the now hydrostatic pressure
517      ! ------------------------------------
518
519      zalfg = 0.5 * grav * rau0
520     
521      zprn(:,:,1) = zalfg * fse3w(:,:,1) * ( 1 + rhd(:,:,1) )       ! Surface value
522
523#if defined key_z_first
524      DO jj = 1, jpj
525         DO ji = 1, jpi
526            DO jk = 2, jpkm1                                        ! Vertical integration from the surface
527               zprn(ji,jj,jk) = zprn(ji,jj,jk-1)   &
528                  &           + zalfg * fse3w(ji,jj,jk) * ( 2. + rhd(ji,jj,jk) + rhd(ji,jj,jk-1) )
529            END DO 
530         END DO 
531      END DO 
532#else
533      DO jk = 2, jpkm1                                              ! Vertical integration from the surface
534         zprn(:,:,jk) = zprn(:,:,jk-1)   &
535            &         + zalfg * fse3w(:,:,jk) * ( 2. + rhd(:,:,jk) + rhd(:,:,jk-1) )
536      END DO 
537#endif
538
539      ! Compute geostrophic balance
540      ! ---------------------------
541#if defined key_z_first
542      DO jj = 2, jpjm1
543         DO ji = 2, jpim1
544            DO jk = 1, jpkm1
545#else
546      DO jk = 1, jpkm1
547         DO jj = 2, jpjm1
548            DO ji = fs_2, fs_jpim1   ! vector opt.
549#endif
550               zmsv = 1. / MAX(  umask(ji-1,jj+1,jk) + umask(ji  ,jj+1,jk)   &
551                               + umask(ji-1,jj  ,jk) + umask(ji  ,jj  ,jk) , 1.  )
552               zphv = ( zprn(ji  ,jj+1,jk) - zprn(ji-1,jj+1,jk) ) * umask(ji-1,jj+1,jk) / e1u(ji-1,jj+1)   &
553                    + ( zprn(ji+1,jj+1,jk) - zprn(ji  ,jj+1,jk) ) * umask(ji  ,jj+1,jk) / e1u(ji  ,jj+1)   &
554                    + ( zprn(ji  ,jj  ,jk) - zprn(ji-1,jj  ,jk) ) * umask(ji-1,jj  ,jk) / e1u(ji-1,jj  )   &
555                    + ( zprn(ji+1,jj  ,jk) - zprn(ji  ,jj  ,jk) ) * umask(ji  ,jj  ,jk) / e1u(ji  ,jj  )
556               zphv = 1. / rau0 * zphv * zmsv * vmask(ji,jj,jk)
557
558               zmsu = 1. / MAX(  vmask(ji+1,jj  ,jk) + vmask(ji  ,jj  ,jk)   &
559                               + vmask(ji+1,jj-1,jk) + vmask(ji  ,jj-1,jk) , 1.  )
560               zphu = ( zprn(ji+1,jj+1,jk) - zprn(ji+1,jj  ,jk) ) * vmask(ji+1,jj  ,jk) / e2v(ji+1,jj  )   &
561                    + ( zprn(ji  ,jj+1,jk) - zprn(ji  ,jj  ,jk) ) * vmask(ji  ,jj  ,jk) / e2v(ji  ,jj  )   &
562                    + ( zprn(ji+1,jj  ,jk) - zprn(ji+1,jj-1,jk) ) * vmask(ji+1,jj-1,jk) / e2v(ji+1,jj-1)   &
563                    + ( zprn(ji  ,jj  ,jk) - zprn(ji  ,jj-1,jk) ) * vmask(ji  ,jj-1,jk) / e2v(ji  ,jj-1)
564               zphu = 1. / rau0 * zphu * zmsu * umask(ji,jj,jk)
565
566               ! Compute the geostrophic velocities
567               un(ji,jj,jk) = -2. * zphu / ( ff(ji,jj) + ff(ji  ,jj-1) )
568               vn(ji,jj,jk) =  2. * zphv / ( ff(ji,jj) + ff(ji-1,jj  ) )
569            END DO
570         END DO
571      END DO
572
573      IF(lwp) WRITE(numout,*) '         we force to zero bottom velocity'
574
575      ! Susbtract the bottom velocity (level jpk-1 for flat bottom case)
576      ! to have a zero bottom velocity
577
578#if defined key_z_first
579      DO jj = 1, jpj
580         DO ji = 1, jpi
581            DO jk = 1, jpkm1
582               un(ji,jj,jk) = ( un(ji,jj,jk) - un(ji,jj,jpkm1) ) * umask(ji,jj,jk)
583               vn(ji,jj,jk) = ( vn(ji,jj,jk) - vn(ji,jj,jpkm1) ) * vmask(ji,jj,jk)
584            END DO
585         END DO
586      END DO
587#else
588      DO jk = 1, jpkm1
589         un(:,:,jk) = ( un(:,:,jk) - un(:,:,jpkm1) ) * umask(:,:,jk)
590         vn(:,:,jk) = ( vn(:,:,jk) - vn(:,:,jpkm1) ) * vmask(:,:,jk)
591      END DO
592#endif
593
594      CALL lbc_lnk( un, 'U', -1. )
595      CALL lbc_lnk( vn, 'V', -1. )
596     
597      ub(:,:,:) = un(:,:,:)
598      vb(:,:,:) = vn(:,:,:)
599     
600      ! WARNING !!!!!
601      ! after initializing u and v, we need to calculate the initial streamfunction bsf.
602      ! Otherwise, only the trend will be computed and the model will blow up (inconsistency).
603      ! to do that, we call dyn_spg with a special trick:
604      ! we fill ua and va with the velocities divided by dt, and the streamfunction will be brought to the
605      ! right value assuming the velocities have been set up in one time step.
606      ! we then set bsfd to zero (first guess for next step is d(psi)/dt = 0.)
607      !  sets up s false trend to calculate the barotropic streamfunction.
608
609      ua(:,:,:) = ub(:,:,:) / rdt
610      va(:,:,:) = vb(:,:,:) / rdt
611
612      ! calls dyn_spg. we assume euler time step, starting from rest.
613      indic = 0
614      CALL dyn_spg( nit000, indic )       ! surface pressure gradient
615
616      ! the new velocity is ua*rdt
617
618      CALL lbc_lnk( ua, 'U', -1. )
619      CALL lbc_lnk( va, 'V', -1. )
620
621      ub(:,:,:) = ua(:,:,:) * rdt
622      vb(:,:,:) = va(:,:,:) * rdt
623      ua(:,:,:) = 0.e0
624      va(:,:,:) = 0.e0
625      un(:,:,:) = ub(:,:,:)
626      vn(:,:,:) = vb(:,:,:)
627       
628      ! Compute the divergence and curl
629
630      CALL div_cur( nit000 )            ! now horizontal divergence and curl
631
632      hdivb(:,:,:) = hdivn(:,:,:)       ! set the before to the now value
633      rotb (:,:,:) = rotn (:,:,:)       ! set the before to the now value
634      !
635      IF( wrk_not_released(3, 1) ) THEN
636         CALL ctl_stop('istate_uvg: failed to release workspace array')
637      ENDIF
638      !
639   END SUBROUTINE istate_uvg
640
641   !!=====================================================================
642END MODULE istate
Note: See TracBrowser for help on using the repository browser.