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 tags/nemo_v1_13_dev2/NEMO/OPA_SRC – NEMO

source: tags/nemo_v1_13_dev2/NEMO/OPA_SRC/istate.F90 @ 2007

Last change on this file since 2007 was 473, checked in by opalod, 18 years ago

nemo_v1_update_060: SM: IOM + 301 levels + CORE + begining of ctl_stop

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 25.3 KB
Line 
1MODULE istate
2   !!======================================================================
3   !!                     ***  MODULE  istate  ***
4   !! Ocean state   :  initial state setting
5   !!=====================================================================
6
7   !!----------------------------------------------------------------------
8   !!   istate_init   : initial state setting
9   !!   istate_tem    : analytical profile for initial Temperature
10   !!   istate_sal    : analytical profile for initial Salinity
11   !!   istate_eel    : initial state setting of EEL R5 configuration
12   !!   istate_gyre   : initial state setting of GYRE configuration
13   !!   istate_uvg    : initial velocity in geostropic balance
14   !!----------------------------------------------------------------------
15   !! * Modules used
16   USE oce             ! ocean dynamics and active tracers
17   USE dom_oce         ! ocean space and time domain
18   USE daymod          !
19   USE ldftra_oce      ! ocean active tracers: lateral physics
20   USE zdf_oce         ! ocean vertical physics
21   USE in_out_manager  ! I/O manager
22   USE phycst          ! physical constants
23   USE wzvmod          ! verctical velocity               (wzv     routine)
24   USE dtatem          ! temperature data                 (dta_tem routine)
25   USE dtasal          ! salinity data                    (dta_sal routine)
26   USE restart         ! ocean restart                   (rst_read routine)
27   USE solisl          ! ???
28
29   IMPLICIT NONE
30   PRIVATE
31
32   !! * Routine accessibility
33   PUBLIC istate_init   ! routine called by step.F90
34
35   !! * Substitutions
36#  include "domzgr_substitute.h90"
37#  include "vectopt_loop_substitute.h90"
38   !!----------------------------------------------------------------------
39   !!   OPA 9.0 , LOCEAN-IPSL (2005)
40   !! $Header$
41   !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt
42   !!----------------------------------------------------------------------
43
44CONTAINS
45
46   SUBROUTINE istate_init
47      !!----------------------------------------------------------------------
48      !!                   ***  ROUTINE istate_init  ***
49      !!
50      !! ** Purpose :   Initialization of the dynamics and tracers.
51      !!
52      !! ** Method  :
53      !!
54      !! History :
55      !!   4.0  !  91-03  ()  Original code
56      !!        !  91-11  (G. Madec)
57      !!   9.0  !  03-09  (G. Madec)  F90: Free form, modules, orthogonality
58      !!----------------------------------------------------------------------
59      USE iom
60      !! * Local declarations
61      !CT INTEGER ::   inum
62      !!----------------------------------------------------------------------
63
64
65      ! Initialization to zero
66      ! ----------------------
67
68      !     before fields       !       now fields        !      after fields       !
69      ;   ub   (:,:,:) = 0.e0   ;   un   (:,:,:) = 0.e0   ;   ua   (:,:,:) = 0.e0
70      ;   vb   (:,:,:) = 0.e0   ;   vn   (:,:,:) = 0.e0   ;   va   (:,:,:) = 0.e0
71      ;                         ;   wn   (:,:,:) = 0.e0   ;
72      ;   rotb (:,:,:) = 0.e0   ;   rotn (:,:,:) = 0.e0   ;
73      ;   hdivb(:,:,:) = 0.e0   ;   hdivn(:,:,:) = 0.e0   ;
74
75      ;   tb   (:,:,:) = 0.e0   ;   tn   (:,:,:) = 0.e0   ;   ta   (:,:,:) = 0.e0
76      ;   sb   (:,:,:) = 0.e0   ;   sn   (:,:,:) = 0.e0   ;   sa   (:,:,:) = 0.e0
77
78      rhd  (:,:,:) = 0.e0
79      rhop (:,:,:) = 0.e0
80      rn2  (:,:,:) = 0.e0 
81
82#if defined key_dynspg_rl
83      ! rigid-lid formulation
84      bsfb(:,:) = 0.e0      ! before barotropic stream-function
85      bsfn(:,:) = 0.e0      ! now    barotropic stream-function
86      bsfd(:,:) = 0.e0      ! barotropic stream-function trend
87#endif
88      ! free surface formulation
89      sshb(:,:) = 0.e0      ! before sea-surface height
90      sshn(:,:) = 0.e0      ! now    sea-surface height
91
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      ELSE
98         !                                    ! Start from rest
99         !                                    ! ---------------
100         neuler = 0                              ! Set time-step indicator at nit000 (euler forward)
101         adatrj = 0._wp
102         IF( cp_cfg == 'eel' ) THEN
103            CALL istate_eel                      ! EEL   configuration : start from pre-defined
104            !                                    !                       velocity and thermohaline fields
105         ELSEIF( cp_cfg == 'gyre' ) THEN         
106            CALL istate_gyre                     ! GYRE  configuration : start from pre-defined temperature
107            !                                    !                       and salinity fields
108         ELSE
109         !                                       ! Other configurations: Initial temperature and salinity fields
110
111         !CT CALL iom_open ('ssh', inum)
112         !CT CALL iom_get( inum, jpdom_local, 'sshb', sshb )     ! free surface formulation (ssh)
113         !CT sshn(:,:) = sshb(:,:)
114         !CT CALL iom_close (inum)
115
116#if defined key_dtatem
117            CALL dta_tem( nit000 )                  ! read 3D temperature data
118            tb(:,:,:) = t_dta(:,:,:)                ! use temperature data read
119            tn(:,:,:) = t_dta(:,:,:)
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(:,:,:)                ! use salinity data read
128            sn(:,:,:) = s_dta(:,:,:)
129#else
130            ! No salinity data
131            IF(lwp)WRITE(numout,*)                  ! analytical salinity profile
132            IF(lwp)WRITE(numout,*)' Salinity initialisation using a constant value'
133            CALL istate_sal
134#endif
135         ENDIF
136
137      ENDIF
138      !                                       ! Vertical velocity
139      !                                       ! -----------------
140      CALL wzv( nit000 )                         ! from horizontal divergence
141
142   END SUBROUTINE istate_init
143
144
145   SUBROUTINE istate_tem
146      !!---------------------------------------------------------------------
147      !!                  ***  ROUTINE istate_tem  ***
148      !!   
149      !! ** Purpose :   Intialization of the temperature field with an
150      !!      analytical profile or a file (i.e. in EEL configuration)
151      !!
152      !! ** Method  :   Use Philander analytic profile of temperature
153      !!
154      !! References :  Philander ???
155      !!
156      !! History :
157      !!   4.0  !  89-12  (P. Andrich)  Original code
158      !!   6.0  !  96-01  (G. Madec)  terrain following coordinates
159      !!   9.0  !  02-09  (G. Madec)  F90: Free form
160      !!----------------------------------------------------------------------
161      !! * Local declarations
162      INTEGER :: ji, jj, jk
163      !!----------------------------------------------------------------------
164
165      IF(lwp) WRITE(numout,*)
166      IF(lwp) WRITE(numout,*) 'istate_tem : initial temperature profile'
167      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
168
169      DO jk = 1, jpk
170         DO jj = 1, jpj
171            DO ji = 1, jpi
172               tn(ji,jj,jk) = (  ( ( 7.5 - 0.*ABS(gphit(ji,jj))/30. )   &
173                  &               *( 1.-TANH((fsdept(ji,jj,jk)-80.)/30.) )   &
174                  &            + 10.*(5000.-fsdept(ji,jj,jk))/5000.)  ) * tmask(ji,jj,jk)
175               tb(ji,jj,jk) = tn(ji,jj,jk)
176          END DO
177        END DO
178      END DO
179
180      IF(lwp) CALL prizre( tn    , jpi   , jpj   , jpk   , jpj/2 ,   &
181         &                 1     , jpi   , 5     , 1     , jpk   ,   &
182         &                 1     , 1.    , numout                  )
183
184   END SUBROUTINE istate_tem
185
186
187   SUBROUTINE istate_sal
188      !!---------------------------------------------------------------------
189      !!                  ***  ROUTINE istate_sal  ***
190      !!
191      !! ** Purpose :   Intialize the salinity field with an analytic profile
192      !!
193      !! ** Method  :   Use to a constant value 35.5
194      !!             
195      !! ** Action  :   Initialize sn and sb
196      !!
197      !! History :
198      !!   4.0  !  89-12  (P. Andrich)  Original code
199      !!   8.5  !  02-09  (G. Madec)  F90: Free form
200      !!----------------------------------------------------------------------
201      !! * Local declarations
202      REAL(wp) ::   zsal = 35.50_wp
203      !!----------------------------------------------------------------------
204
205      IF(lwp) WRITE(numout,*)
206      IF(lwp) WRITE(numout,*) 'istate_sal : initial salinity : ', zsal
207      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
208
209      sn(:,:,:) = zsal * tmask(:,:,:)
210      sb(:,:,:) = sn(:,:,:)
211     
212   END SUBROUTINE istate_sal
213
214
215   SUBROUTINE istate_eel
216      !!----------------------------------------------------------------------
217      !!                   ***  ROUTINE istate_eel  ***
218      !!
219      !! ** Purpose :   Initialization of the dynamics and tracers for EEL R5
220      !!      configuration (channel with or without a topographic bump)
221      !!
222      !! ** Method  : - set temprature field
223      !!              - set salinity field
224      !!              - set velocity field including horizontal divergence
225      !!                and relative vorticity fields
226      !!
227      !! History :
228      !!   8.0  !  01-09  (M. Levy, M. Ben Jelloul)  read file for EEL 2 & 6
229      !!   9.0  !  03-09  (G. Madec, C. Talandier)  EEL 5
230      !!   9.0  !  05-11  (V. Garnier) Surface pressure gradient organization
231      !!----------------------------------------------------------------------
232      !! * Modules used
233      USE eosbn2     ! eq. of state, Brunt Vaisala frequency (eos     routine)
234      USE divcur     ! hor. divergence & rel. vorticity      (div_cur routine)
235      USE iom
236 
237      !! * Local declarations
238      INTEGER  ::   inum              ! temporary logical unit
239      INTEGER  ::   ji, jj, jk        ! dummy loop indices
240      REAL(wp) ::   &
241         zh1, zh2, zslope, zcst       ! temporary scalars
242      REAL(wp) ::   &
243         zt1  = 12._wp,            &  ! surface temperature value (EEL R5)
244         zt2  =  2._wp,            &  ! bottom  temperature value (EEL R5)
245         zsal = 35.5_wp                ! constant salinity (EEL R2, R5 and R6)
246# if ! defined key_dynspg_rl
247      REAL(wp), DIMENSION(jpiglo,jpjglo) ::   &
248         zssh                         ! initial ssh over the global domain
249# endif
250      !!----------------------------------------------------------------------
251
252      SELECT CASE ( jp_cfg ) 
253         !                                              ! ====================
254         CASE ( 5 )                                     ! EEL R5 configuration
255            !                                           ! ====================
256
257            ! set temperature field with a linear profile
258            ! -------------------------------------------
259            IF(lwp) WRITE(numout,*)
260            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: linear temperature profile'
261            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
262
263            zh1 = gdept_0(  1  )
264            zh2 = gdept_0(jpkm1)
265
266            zslope = ( zt1 - zt2 ) / ( zh1 - zh2 )
267            zcst   = ( zt1 * ( zh1 - zh2) - ( zt1 - zt2 ) * zh1 ) / ( zh1 - zh2 )
268
269            DO jk = 1, jpk
270               tn(:,:,jk) = ( zslope * fsdept(:,:,jk) + zcst  ) * tmask(:,:,jk)
271               tb(:,:,jk) = tn(:,:,jk)
272            END DO
273
274            IF(lwp) CALL prizre( tn    , jpi   , jpj   , jpk   , jpj/2 ,   &
275               &                 1     , jpi   , 5     , 1     , jpk   ,   &
276               &                 1     , 1.    , numout                  )
277
278            ! set salinity field to a constant value
279            ! --------------------------------------
280            IF(lwp) WRITE(numout,*)
281            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal
282            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
283
284            sn(:,:,:) = zsal * tmask(:,:,:)
285            sb(:,:,:) = sn(:,:,:)
286     
287
288# if ! defined key_dynspg_rl
289            ! set the dynamics: U,V, hdiv, rot (and ssh if necessary)
290            ! ----------------
291            ! Start EEL5 configuration with barotropic geostrophic velocities
292            ! according the sshb and sshn SSH imposed.
293            ub(:,:,:) = 0.1 * umask(:,:,:)
294            un(:,:,:) = ub(:,:,:)
295
296            DO jj = 1, jpjglo
297               zssh(:,jj) = ( .22 - ( FLOAT(jj-3) * (0.44) ) / 99. )
298            END DO
299            DO jj = 1, nlcj
300               DO ji = 1, nlci
301                  sshb(ji,jj) = zssh( mig(ji) , mjg(jj) ) * tmask(ji,jj,1)
302               END DO
303            END DO
304            sshb(nlci+1:jpi,      :   ) = 0.e0      ! set to zero extra mpp columns
305            sshb(      :   ,nlcj+1:jpj) = 0.e0      ! set to zero extra mpp rows
306
307            sshn(:,:) = sshb(:,:)                   ! set now ssh to the before value
308
309            ! horizontal divergence and relative vorticity (curl)
310            CALL div_cur( nit000 )
311
312            ! N.B. the vertical velocity will be computed from the horizontal divergence field
313            ! in istate by a call to wzv routine
314# endif
315
316
317            !                                     ! ==========================
318         CASE ( 2 , 6 )                           ! EEL R2 or R6 configuration
319            !                                     ! ==========================
320 
321            ! set temperature field with a NetCDF file
322            ! ----------------------------------------
323            IF(lwp) WRITE(numout,*)
324            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R2 or R6: read initial temperature in a NetCDF file'
325            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
326
327            CALL iom_open ( 'eel.initemp', inum )
328            CALL iom_get ( inum, jpdom_data, 'initemp', tb ) ! read before temprature (tb)
329            CALL iom_close( inum )
330     
331            tn(:,:,:) = tb(:,:,:)                            ! set nox temperature to tb
332
333            IF(lwp) CALL prizre( tn    , jpi   , jpj   , jpk   , jpj/2 ,   &
334               &                 1     , jpi   , 5     , 1     , jpk   ,   &
335               &                 1     , 1.    , numout                  )
336
337
338            ! set salinity field to a constant value
339            ! --------------------------------------
340            IF(lwp) WRITE(numout,*)
341            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal
342            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
343 
344            sn(:,:,:) = zsal * tmask(:,:,:)
345            sb(:,:,:) = sn(:,:,:)
346
347
348            IF( lk_isl ) THEN
349               ! Horizontal velocity : start from geostrophy (EEL config)
350               CALL eos( tn, sn, rhd )     ! now in situ density
351               CALL istate_uvg             ! compute geostrophic velocity
352
353               ! N.B. the vertical velocity will be computed from the horizontal divergence field
354               ! in istate by a call to wzv routine
355            ENDIF
356            !                                    ! ===========================
357         CASE DEFAULT                            ! NONE existing configuration
358            !                                    ! ===========================
359            WRITE(ctmp1,*) 'EEL with a ', jp_cfg,' km resolution is not coded'
360            CALL ctl_stop( ctmp1 )
361
362      END SELECT
363
364   END SUBROUTINE istate_eel
365
366
367   SUBROUTINE istate_gyre
368      !!----------------------------------------------------------------------
369      !!                   ***  ROUTINE istate_gyre  ***
370      !!
371      !! ** Purpose :   Initialization of the dynamics and tracers for GYRE
372      !!      configuration (double gyre with rotated domain)
373      !!
374      !! ** Method  : - set temprature field
375      !!              - set salinity field
376      !!
377      !! ** History :                                     
378      !!      9.0  !  04-05  (A. Koch-Larrouy)  Original code
379      !!----------------------------------------------------------------------
380      !! * Modules used
381      USE iom
382
383      !! * Local variables
384      INTEGER  ::   inum              ! temporary logical unit
385      INTEGER, PARAMETER ::   &
386         ntsinit = 0         ! (0/1) (analytical/input data files) T&S initialization
387
388      INTEGER :: ji, jj, jk  ! dummy loop indices
389      !!----------------------------------------------------------------------
390
391      SELECT CASE ( ntsinit)
392
393      CASE ( 0 )                  ! analytical T/S profil deduced from LEVITUS
394         IF(lwp) WRITE(numout,*)
395         IF(lwp) WRITE(numout,*) 'istate_gyre : initial analytical T and S profil deduced from LEVITUS '
396         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
397
398         DO jk = 1, jpk
399            DO jj = 1, jpj
400               DO ji = 1, jpi
401                  tn(ji,jj,jk) = (  16. - 12. * TANH( (fsdept(ji,jj,jk) - 400) / 700 )         )   &
402                       &           * (-TANH( (500-fsdept(ji,jj,jk)) / 150 ) + 1) / 2               &
403                       &       + (      15. * ( 1. - TANH( (fsdept(ji,jj,jk)-50.) / 1500.) )       &
404                       &                - 1.4 * TANH((fsdept(ji,jj,jk)-100.) / 100.)               &   
405                       &                + 7.  * (1500. - fsdept(ji,jj,jk)) / 1500.             )   & 
406                       &           * (-TANH( (fsdept(ji,jj,jk) - 500) / 150) + 1) / 2
407                  tn(ji,jj,jk) = tn(ji,jj,jk) * tmask(ji,jj,jk)
408                  tb(ji,jj,jk) = tn(ji,jj,jk)
409
410                  sn(ji,jj,jk) =  (  36.25 - 1.13 * TANH( (fsdept(ji,jj,jk) - 305) / 460 )  )  &
411                     &              * (-TANH((500 - fsdept(ji,jj,jk)) / 150) + 1) / 2          &
412                     &          + (  35.55 + 1.25 * (5000. - fsdept(ji,jj,jk)) / 5000.         &
413                     &                - 1.62 * TANH( (fsdept(ji,jj,jk) - 60.  ) / 650. )       &
414                     &                + 0.2  * TANH( (fsdept(ji,jj,jk) - 35.  ) / 100. )       &
415                     &                + 0.2  * TANH( (fsdept(ji,jj,jk) - 1000.) / 5000.)    )  &
416                     &              * (-TANH((fsdept(ji,jj,jk) - 500) / 150) + 1) / 2 
417                  sn(ji,jj,jk) = sn(ji,jj,jk) * tmask(ji,jj,jk)
418                  sb(ji,jj,jk) = sn(ji,jj,jk)
419               END DO
420            END DO
421         END DO
422
423      CASE ( 1 )                  ! T/S data fields read in dta_tem.nc/data_sal.nc files
424         IF(lwp) WRITE(numout,*)
425         IF(lwp) WRITE(numout,*) 'istate_gyre : initial T and S read from dta_tem.nc/data_sal.nc files'
426         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
427         IF(lwp) WRITE(numout,*) '              NetCDF FORMAT'
428
429         ! Read temperature field
430         ! ----------------------
431         CALL iom_open ( 'data_tem', inum )
432         CALL iom_get ( inum, jpdom_data, 'votemper', tn ) 
433         CALL iom_close( inum )
434
435         tn(:,:,:) = tn(:,:,:) * tmask(:,:,:) 
436         tb(:,:,:) = tn(:,:,:)
437
438         ! Read salinity field
439         ! -------------------
440         CALL iom_open ( 'data_sal', inum )
441         CALL iom_get ( inum, jpdom_data, 'vosaline', sn ) 
442         CALL iom_close( inum )
443
444         sn(:,:,:)  = sn(:,:,:) * tmask(:,:,:) 
445         sb(:,:,:)  = sn(:,:,:)
446
447      END SELECT
448
449      IF(lwp) THEN
450         WRITE(numout,*)
451         WRITE(numout,*) '              Initial temperature and salinity profiles:'
452         WRITE(numout, "(9x,' level   gdept_0   temperature   salinity   ')" )
453         WRITE(numout, "(10x, i4, 3f10.2)" ) ( jk, gdept_0(jk), tn(2,2,jk), sn(2,2,jk), jk = 1, jpk )
454      ENDIF
455
456     
457   END SUBROUTINE istate_gyre
458
459
460
461   SUBROUTINE istate_uvg
462      !!----------------------------------------------------------------------
463      !!                  ***  ROUTINE istate_uvg  ***
464      !!
465      !! ** Purpose :   Compute the geostrophic velocities from (tn,sn) fields
466      !!
467      !! ** Method  :   Using the hydrostatic hypothesis the now hydrostatic
468      !!      pressure is computed by integrating the in-situ density from the
469      !!      surface to the bottom.
470      !!                 p=integral [ rau*g dz ]
471      !!
472      !! History :
473      !!   8.1  !  01-09  (M. Levy, M. Ben Jelloul)  Original code
474      !!   8.5  !  02-09  (G. Madec)  F90: Free form
475      !!   9.0  !  05-11  (V. Garnier) Surface pressure gradient organization
476      !!----------------------------------------------------------------------
477      !! * Modules used
478      USE eosbn2          ! eq. of state, Brunt Vaisala frequency (eos     routine)
479      USE dynspg          ! surface pressure gradient             (dyn_spg routine)
480      USE divcur          ! hor. divergence & rel. vorticity      (div_cur routine)
481      USE lbclnk          ! ocean lateral boundary condition (or mpp link)
482
483      !! * Local declarations
484      INTEGER ::   ji, jj, jk        ! dummy loop indices
485      INTEGER ::   indic             ! ???
486      REAL(wp) ::   &
487         zmsv, zphv, zmsu, zphu,  &  ! temporary scalars
488         zalfg
489      REAL(wp), DIMENSION (jpi,jpj,jpk) ::   &
490         zprn                        ! workspace
491      !!----------------------------------------------------------------------
492
493      IF(lwp) WRITE(numout,*) 
494      IF(lwp) WRITE(numout,*) 'istate_uvg : Start from Geostrophy'
495      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
496
497      ! Compute the now hydrostatic pressure
498      ! ------------------------------------
499
500      zalfg = 0.5 * grav * rau0
501      ! Surface value
502      zprn(:,:,1) = zalfg * fse3w(:,:,1) * ( 1 + rhd(:,:,1) )
503
504      ! Vertical integration from the surface
505      DO jk = 2, jpkm1
506         zprn(:,:,jk) = zprn(:,:,jk-1)   &
507            &         + zalfg * fse3w(:,:,jk) * ( 2. + rhd(:,:,jk) + rhd(:,:,jk-1) )
508      END DO 
509
510      ! Compute geostrophic balance
511      ! ---------------------------
512
513      DO jk = 1, jpkm1
514         DO jj = 2, jpjm1
515            DO ji = fs_2, fs_jpim1   ! vertor opt.
516               zmsv = 1. / MAX(  umask(ji-1,jj+1,jk) + umask(ji  ,jj+1,jk)   &
517                               + umask(ji-1,jj  ,jk) + umask(ji  ,jj  ,jk) , 1.  )
518               zphv = ( zprn(ji  ,jj+1,jk) - zprn(ji-1,jj+1,jk) ) * umask(ji-1,jj+1,jk) / e1u(ji-1,jj+1)   &
519                    + ( zprn(ji+1,jj+1,jk) - zprn(ji  ,jj+1,jk) ) * umask(ji  ,jj+1,jk) / e1u(ji  ,jj+1)   &
520                    + ( zprn(ji  ,jj  ,jk) - zprn(ji-1,jj  ,jk) ) * umask(ji-1,jj  ,jk) / e1u(ji-1,jj  )   &
521                    + ( zprn(ji+1,jj  ,jk) - zprn(ji  ,jj  ,jk) ) * umask(ji  ,jj  ,jk) / e1u(ji  ,jj  )
522               zphv = 1. / rau0 * zphv * zmsv * vmask(ji,jj,jk)
523
524               zmsu = 1. / MAX(  vmask(ji+1,jj  ,jk) + vmask(ji  ,jj  ,jk)   &
525                               + vmask(ji+1,jj-1,jk) + vmask(ji  ,jj-1,jk) , 1.  )
526               zphu = ( zprn(ji+1,jj+1,jk) - zprn(ji+1,jj  ,jk) ) * vmask(ji+1,jj  ,jk) / e2v(ji+1,jj  )   &
527                    + ( zprn(ji  ,jj+1,jk) - zprn(ji  ,jj  ,jk) ) * vmask(ji  ,jj  ,jk) / e2v(ji  ,jj  )   &
528                    + ( zprn(ji+1,jj  ,jk) - zprn(ji+1,jj-1,jk) ) * vmask(ji+1,jj-1,jk) / e2v(ji+1,jj-1)   &
529                    + ( zprn(ji  ,jj  ,jk) - zprn(ji  ,jj-1,jk) ) * vmask(ji  ,jj-1,jk) / e2v(ji  ,jj-1)
530               zphu = 1. / rau0 * zphu * zmsu * umask(ji,jj,jk)
531
532               ! Compute the geostrophic velocities
533               un(ji,jj,jk) = -2. * zphu / ( ff(ji,jj) + ff(ji  ,jj-1) )
534               vn(ji,jj,jk) =  2. * zphv / ( ff(ji,jj) + ff(ji-1,jj  ) )
535            END DO
536         END DO
537      END DO
538
539      IF(lwp) WRITE(numout,*) '         we force to zero bottom velocity'
540
541      ! Susbtract the bottom velocity (level jpk-1 for flat bottom case)
542      ! to have a zero bottom velocity
543
544      DO jk = 1, jpkm1
545         un(:,:,jk) = ( un(:,:,jk) - un(:,:,jpkm1) ) * umask(:,:,jk)
546         vn(:,:,jk) = ( vn(:,:,jk) - vn(:,:,jpkm1) ) * vmask(:,:,jk)
547      END DO
548
549      CALL lbc_lnk( un, 'U', -1. )
550      CALL lbc_lnk( vn, 'V', -1. )
551     
552      ub(:,:,:) = un(:,:,:)
553      vb(:,:,:) = vn(:,:,:)
554     
555      ! WARNING !!!!!
556      ! after initializing u and v, we need to calculate the initial streamfunction bsf.
557      ! Otherwise, only the trend will be computed and the model will blow up (inconsistency).
558     
559      ! to do that, we call dyn_spg with a special trick:
560      ! we fill ua and va with the velocities divided by dt,
561      ! and the streamfunction will be brought to the right
562      ! value assuming the velocities have been set up in
563      ! one time step.
564      ! we then set bsfd to zero (first guess for next step
565      ! is d(psi)/dt = 0.)
566
567      !  sets up s false trend to calculate the barotropic
568      !  streamfunction.
569
570      ua(:,:,:) = ub(:,:,:) / rdt
571      va(:,:,:) = vb(:,:,:) / rdt
572
573      ! calls dyn_spg. we assume euler time step, starting from rest.
574      indic = 0
575      CALL dyn_spg( nit000, indic )       ! surface pressure gradient
576
577      ! the new velocity is ua*rdt
578
579      CALL lbc_lnk( ua, 'U', -1. )
580      CALL lbc_lnk( va, 'V', -1. )
581
582      ub(:,:,:) = ua(:,:,:) * rdt
583      vb(:,:,:) = va(:,:,:) * rdt
584      ua(:,:,:) = 0.e0
585      va(:,:,:) = 0.e0
586      un(:,:,:) = ub(:,:,:)
587      vn(:,:,:) = vb(:,:,:)
588       
589#if defined key_dynspg_rl
590      IF( lk_isl )   bsfb(:,:) = bsfn(:,:)          ! Put bsfb to zero
591#endif
592
593      ! Compute the divergence and curl
594
595      CALL div_cur( nit000 )            ! now horizontal divergence and curl
596
597      hdivb(:,:,:) = hdivn(:,:,:)       ! set the before to the now value
598      rotb (:,:,:) = rotn (:,:,:)       ! set the before to the now value
599
600   END SUBROUTINE istate_uvg
601
602   !!=====================================================================
603END MODULE istate
Note: See TracBrowser for help on using the repository browser.