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

source: trunk/NEMO/OPA_SRC/istate.F90 @ 507

Last change on this file since 507 was 479, checked in by opalod, 18 years ago

nemo_v1_bugfix_048 : CT : correction of the ssh initial state computation

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 25.9 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      INTEGER  ::   ijloc
241      REAL(wp) ::   &
242         zh1, zh2, zslope, zcst,   &  ! temporary scalars
243         zfcor
244      REAL(wp) ::   &
245         zt1  = 12._wp,            &  ! surface temperature value (EEL R5)
246         zt2  =  2._wp,            &  ! bottom  temperature value (EEL R5)
247         zsal = 35.5_wp,           &  ! constant salinity (EEL R2, R5 and R6)
248         zueel = 0.1_wp               ! constant uniform zonal velocity (EEL R5)
249# if ! defined key_dynspg_rl
250      REAL(wp), DIMENSION(jpiglo,jpjglo) ::   &
251         zssh                         ! initial ssh over the global domain
252# endif
253      !!----------------------------------------------------------------------
254
255      SELECT CASE ( jp_cfg ) 
256         !                                              ! ====================
257         CASE ( 5 )                                     ! EEL R5 configuration
258            !                                           ! ====================
259
260            ! set temperature field with a linear profile
261            ! -------------------------------------------
262            IF(lwp) WRITE(numout,*)
263            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: linear temperature profile'
264            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
265
266            zh1 = gdept_0(  1  )
267            zh2 = gdept_0(jpkm1)
268
269            zslope = ( zt1 - zt2 ) / ( zh1 - zh2 )
270            zcst   = ( zt1 * ( zh1 - zh2) - ( zt1 - zt2 ) * zh1 ) / ( zh1 - zh2 )
271
272            DO jk = 1, jpk
273               tn(:,:,jk) = ( zslope * fsdept(:,:,jk) + zcst  ) * tmask(:,:,jk)
274               tb(:,:,jk) = tn(:,:,jk)
275            END DO
276
277            IF(lwp) CALL prizre( tn    , jpi   , jpj   , jpk   , jpj/2 ,   &
278               &                 1     , jpi   , 5     , 1     , jpk   ,   &
279               &                 1     , 1.    , numout                  )
280
281            ! set salinity field to a constant value
282            ! --------------------------------------
283            IF(lwp) WRITE(numout,*)
284            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal
285            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
286
287            sn(:,:,:) = zsal * tmask(:,:,:)
288            sb(:,:,:) = sn(:,:,:)
289     
290
291# if ! defined key_dynspg_rl
292            ! set the dynamics: U,V, hdiv, rot (and ssh if necessary)
293            ! ----------------
294            ! Start EEL5 configuration with barotropic geostrophic velocities
295            ! according the sshb and sshn SSH imposed.
296            ! we assume a uniform grid (hence the use of e1t(1,1) for delta_y)
297            ! we use the Coriolis frequency at mid-channel.   
298   
299            ub(:,:,:) = zueel * umask(:,:,:)
300            un(:,:,:) = ub(:,:,:)
301            ijloc = mj0(INT(jpjglo-1)/2)
302            zfcor = ff(1,ijloc)
303
304            DO jj = 1, jpjglo
305               zssh(:,jj) = - (FLOAT(jj)- FLOAT(jpjglo-1)/2.)*zueel*e1t(1,1)*zfcor/grav 
306            END DO
307
308            IF(lwp) THEN
309               WRITE(numout,*) ' Uniform zonal velocity for EEL R5:',zueel
310               WRITE(numout,*) ' Geostrophic SSH profile as a function of y:'
311               WRITE(numout,'(12(1x,f6.2))') zssh(1,:)
312            ENDIF
313
314            DO jj = 1, nlcj
315               DO ji = 1, nlci
316                  sshb(ji,jj) = zssh( mig(ji) , mjg(jj) ) * tmask(ji,jj,1)
317               END DO
318            END DO
319            sshb(nlci+1:jpi,      :   ) = 0.e0      ! set to zero extra mpp columns
320            sshb(      :   ,nlcj+1:jpj) = 0.e0      ! set to zero extra mpp rows
321
322            sshn(:,:) = sshb(:,:)                   ! set now ssh to the before value
323
324            ! horizontal divergence and relative vorticity (curl)
325            CALL div_cur( nit000 )
326
327            ! N.B. the vertical velocity will be computed from the horizontal divergence field
328            ! in istate by a call to wzv routine
329# endif
330
331
332            !                                     ! ==========================
333         CASE ( 2 , 6 )                           ! EEL R2 or R6 configuration
334            !                                     ! ==========================
335 
336            ! set temperature field with a NetCDF file
337            ! ----------------------------------------
338            IF(lwp) WRITE(numout,*)
339            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R2 or R6: read initial temperature in a NetCDF file'
340            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
341
342            CALL iom_open ( 'eel.initemp', inum )
343            CALL iom_get ( inum, jpdom_data, 'initemp', tb ) ! read before temprature (tb)
344            CALL iom_close( inum )
345     
346            tn(:,:,:) = tb(:,:,:)                            ! set nox temperature to tb
347
348            IF(lwp) CALL prizre( tn    , jpi   , jpj   , jpk   , jpj/2 ,   &
349               &                 1     , jpi   , 5     , 1     , jpk   ,   &
350               &                 1     , 1.    , numout                  )
351
352
353            ! set salinity field to a constant value
354            ! --------------------------------------
355            IF(lwp) WRITE(numout,*)
356            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal
357            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
358 
359            sn(:,:,:) = zsal * tmask(:,:,:)
360            sb(:,:,:) = sn(:,:,:)
361
362
363            IF( lk_isl ) THEN
364               ! Horizontal velocity : start from geostrophy (EEL config)
365               CALL eos( tn, sn, rhd )     ! now in situ density
366               CALL istate_uvg             ! compute geostrophic velocity
367
368               ! N.B. the vertical velocity will be computed from the horizontal divergence field
369               ! in istate by a call to wzv routine
370            ENDIF
371            !                                    ! ===========================
372         CASE DEFAULT                            ! NONE existing configuration
373            !                                    ! ===========================
374            WRITE(ctmp1,*) 'EEL with a ', jp_cfg,' km resolution is not coded'
375            CALL ctl_stop( ctmp1 )
376
377      END SELECT
378
379   END SUBROUTINE istate_eel
380
381
382   SUBROUTINE istate_gyre
383      !!----------------------------------------------------------------------
384      !!                   ***  ROUTINE istate_gyre  ***
385      !!
386      !! ** Purpose :   Initialization of the dynamics and tracers for GYRE
387      !!      configuration (double gyre with rotated domain)
388      !!
389      !! ** Method  : - set temprature field
390      !!              - set salinity field
391      !!
392      !! ** History :                                     
393      !!      9.0  !  04-05  (A. Koch-Larrouy)  Original code
394      !!----------------------------------------------------------------------
395      !! * Modules used
396      USE iom
397
398      !! * Local variables
399      INTEGER  ::   inum              ! temporary logical unit
400      INTEGER, PARAMETER ::   &
401         ntsinit = 0         ! (0/1) (analytical/input data files) T&S initialization
402
403      INTEGER :: ji, jj, jk  ! dummy loop indices
404      !!----------------------------------------------------------------------
405
406      SELECT CASE ( ntsinit)
407
408      CASE ( 0 )                  ! analytical T/S profil deduced from LEVITUS
409         IF(lwp) WRITE(numout,*)
410         IF(lwp) WRITE(numout,*) 'istate_gyre : initial analytical T and S profil deduced from LEVITUS '
411         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
412
413         DO jk = 1, jpk
414            DO jj = 1, jpj
415               DO ji = 1, jpi
416                  tn(ji,jj,jk) = (  16. - 12. * TANH( (fsdept(ji,jj,jk) - 400) / 700 )         )   &
417                       &           * (-TANH( (500-fsdept(ji,jj,jk)) / 150 ) + 1) / 2               &
418                       &       + (      15. * ( 1. - TANH( (fsdept(ji,jj,jk)-50.) / 1500.) )       &
419                       &                - 1.4 * TANH((fsdept(ji,jj,jk)-100.) / 100.)               &   
420                       &                + 7.  * (1500. - fsdept(ji,jj,jk)) / 1500.             )   & 
421                       &           * (-TANH( (fsdept(ji,jj,jk) - 500) / 150) + 1) / 2
422                  tn(ji,jj,jk) = tn(ji,jj,jk) * tmask(ji,jj,jk)
423                  tb(ji,jj,jk) = tn(ji,jj,jk)
424
425                  sn(ji,jj,jk) =  (  36.25 - 1.13 * TANH( (fsdept(ji,jj,jk) - 305) / 460 )  )  &
426                     &              * (-TANH((500 - fsdept(ji,jj,jk)) / 150) + 1) / 2          &
427                     &          + (  35.55 + 1.25 * (5000. - fsdept(ji,jj,jk)) / 5000.         &
428                     &                - 1.62 * TANH( (fsdept(ji,jj,jk) - 60.  ) / 650. )       &
429                     &                + 0.2  * TANH( (fsdept(ji,jj,jk) - 35.  ) / 100. )       &
430                     &                + 0.2  * TANH( (fsdept(ji,jj,jk) - 1000.) / 5000.)    )  &
431                     &              * (-TANH((fsdept(ji,jj,jk) - 500) / 150) + 1) / 2 
432                  sn(ji,jj,jk) = sn(ji,jj,jk) * tmask(ji,jj,jk)
433                  sb(ji,jj,jk) = sn(ji,jj,jk)
434               END DO
435            END DO
436         END DO
437
438      CASE ( 1 )                  ! T/S data fields read in dta_tem.nc/data_sal.nc files
439         IF(lwp) WRITE(numout,*)
440         IF(lwp) WRITE(numout,*) 'istate_gyre : initial T and S read from dta_tem.nc/data_sal.nc files'
441         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
442         IF(lwp) WRITE(numout,*) '              NetCDF FORMAT'
443
444         ! Read temperature field
445         ! ----------------------
446         CALL iom_open ( 'data_tem', inum )
447         CALL iom_get ( inum, jpdom_data, 'votemper', tn ) 
448         CALL iom_close( inum )
449
450         tn(:,:,:) = tn(:,:,:) * tmask(:,:,:) 
451         tb(:,:,:) = tn(:,:,:)
452
453         ! Read salinity field
454         ! -------------------
455         CALL iom_open ( 'data_sal', inum )
456         CALL iom_get ( inum, jpdom_data, 'vosaline', sn ) 
457         CALL iom_close( inum )
458
459         sn(:,:,:)  = sn(:,:,:) * tmask(:,:,:) 
460         sb(:,:,:)  = sn(:,:,:)
461
462      END SELECT
463
464      IF(lwp) THEN
465         WRITE(numout,*)
466         WRITE(numout,*) '              Initial temperature and salinity profiles:'
467         WRITE(numout, "(9x,' level   gdept_0   temperature   salinity   ')" )
468         WRITE(numout, "(10x, i4, 3f10.2)" ) ( jk, gdept_0(jk), tn(2,2,jk), sn(2,2,jk), jk = 1, jpk )
469      ENDIF
470
471     
472   END SUBROUTINE istate_gyre
473
474
475
476   SUBROUTINE istate_uvg
477      !!----------------------------------------------------------------------
478      !!                  ***  ROUTINE istate_uvg  ***
479      !!
480      !! ** Purpose :   Compute the geostrophic velocities from (tn,sn) fields
481      !!
482      !! ** Method  :   Using the hydrostatic hypothesis the now hydrostatic
483      !!      pressure is computed by integrating the in-situ density from the
484      !!      surface to the bottom.
485      !!                 p=integral [ rau*g dz ]
486      !!
487      !! History :
488      !!   8.1  !  01-09  (M. Levy, M. Ben Jelloul)  Original code
489      !!   8.5  !  02-09  (G. Madec)  F90: Free form
490      !!   9.0  !  05-11  (V. Garnier) Surface pressure gradient organization
491      !!----------------------------------------------------------------------
492      !! * Modules used
493      USE eosbn2          ! eq. of state, Brunt Vaisala frequency (eos     routine)
494      USE dynspg          ! surface pressure gradient             (dyn_spg routine)
495      USE divcur          ! hor. divergence & rel. vorticity      (div_cur routine)
496      USE lbclnk          ! ocean lateral boundary condition (or mpp link)
497
498      !! * Local declarations
499      INTEGER ::   ji, jj, jk        ! dummy loop indices
500      INTEGER ::   indic             ! ???
501      REAL(wp) ::   &
502         zmsv, zphv, zmsu, zphu,  &  ! temporary scalars
503         zalfg
504      REAL(wp), DIMENSION (jpi,jpj,jpk) ::   &
505         zprn                        ! workspace
506      !!----------------------------------------------------------------------
507
508      IF(lwp) WRITE(numout,*) 
509      IF(lwp) WRITE(numout,*) 'istate_uvg : Start from Geostrophy'
510      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
511
512      ! Compute the now hydrostatic pressure
513      ! ------------------------------------
514
515      zalfg = 0.5 * grav * rau0
516      ! Surface value
517      zprn(:,:,1) = zalfg * fse3w(:,:,1) * ( 1 + rhd(:,:,1) )
518
519      ! Vertical integration from the surface
520      DO jk = 2, jpkm1
521         zprn(:,:,jk) = zprn(:,:,jk-1)   &
522            &         + zalfg * fse3w(:,:,jk) * ( 2. + rhd(:,:,jk) + rhd(:,:,jk-1) )
523      END DO 
524
525      ! Compute geostrophic balance
526      ! ---------------------------
527
528      DO jk = 1, jpkm1
529         DO jj = 2, jpjm1
530            DO ji = fs_2, fs_jpim1   ! vertor opt.
531               zmsv = 1. / MAX(  umask(ji-1,jj+1,jk) + umask(ji  ,jj+1,jk)   &
532                               + umask(ji-1,jj  ,jk) + umask(ji  ,jj  ,jk) , 1.  )
533               zphv = ( zprn(ji  ,jj+1,jk) - zprn(ji-1,jj+1,jk) ) * umask(ji-1,jj+1,jk) / e1u(ji-1,jj+1)   &
534                    + ( zprn(ji+1,jj+1,jk) - zprn(ji  ,jj+1,jk) ) * umask(ji  ,jj+1,jk) / e1u(ji  ,jj+1)   &
535                    + ( zprn(ji  ,jj  ,jk) - zprn(ji-1,jj  ,jk) ) * umask(ji-1,jj  ,jk) / e1u(ji-1,jj  )   &
536                    + ( zprn(ji+1,jj  ,jk) - zprn(ji  ,jj  ,jk) ) * umask(ji  ,jj  ,jk) / e1u(ji  ,jj  )
537               zphv = 1. / rau0 * zphv * zmsv * vmask(ji,jj,jk)
538
539               zmsu = 1. / MAX(  vmask(ji+1,jj  ,jk) + vmask(ji  ,jj  ,jk)   &
540                               + vmask(ji+1,jj-1,jk) + vmask(ji  ,jj-1,jk) , 1.  )
541               zphu = ( zprn(ji+1,jj+1,jk) - zprn(ji+1,jj  ,jk) ) * vmask(ji+1,jj  ,jk) / e2v(ji+1,jj  )   &
542                    + ( zprn(ji  ,jj+1,jk) - zprn(ji  ,jj  ,jk) ) * vmask(ji  ,jj  ,jk) / e2v(ji  ,jj  )   &
543                    + ( zprn(ji+1,jj  ,jk) - zprn(ji+1,jj-1,jk) ) * vmask(ji+1,jj-1,jk) / e2v(ji+1,jj-1)   &
544                    + ( zprn(ji  ,jj  ,jk) - zprn(ji  ,jj-1,jk) ) * vmask(ji  ,jj-1,jk) / e2v(ji  ,jj-1)
545               zphu = 1. / rau0 * zphu * zmsu * umask(ji,jj,jk)
546
547               ! Compute the geostrophic velocities
548               un(ji,jj,jk) = -2. * zphu / ( ff(ji,jj) + ff(ji  ,jj-1) )
549               vn(ji,jj,jk) =  2. * zphv / ( ff(ji,jj) + ff(ji-1,jj  ) )
550            END DO
551         END DO
552      END DO
553
554      IF(lwp) WRITE(numout,*) '         we force to zero bottom velocity'
555
556      ! Susbtract the bottom velocity (level jpk-1 for flat bottom case)
557      ! to have a zero bottom velocity
558
559      DO jk = 1, jpkm1
560         un(:,:,jk) = ( un(:,:,jk) - un(:,:,jpkm1) ) * umask(:,:,jk)
561         vn(:,:,jk) = ( vn(:,:,jk) - vn(:,:,jpkm1) ) * vmask(:,:,jk)
562      END DO
563
564      CALL lbc_lnk( un, 'U', -1. )
565      CALL lbc_lnk( vn, 'V', -1. )
566     
567      ub(:,:,:) = un(:,:,:)
568      vb(:,:,:) = vn(:,:,:)
569     
570      ! WARNING !!!!!
571      ! after initializing u and v, we need to calculate the initial streamfunction bsf.
572      ! Otherwise, only the trend will be computed and the model will blow up (inconsistency).
573     
574      ! to do that, we call dyn_spg with a special trick:
575      ! we fill ua and va with the velocities divided by dt,
576      ! and the streamfunction will be brought to the right
577      ! value assuming the velocities have been set up in
578      ! one time step.
579      ! we then set bsfd to zero (first guess for next step
580      ! is d(psi)/dt = 0.)
581
582      !  sets up s false trend to calculate the barotropic
583      !  streamfunction.
584
585      ua(:,:,:) = ub(:,:,:) / rdt
586      va(:,:,:) = vb(:,:,:) / rdt
587
588      ! calls dyn_spg. we assume euler time step, starting from rest.
589      indic = 0
590      CALL dyn_spg( nit000, indic )       ! surface pressure gradient
591
592      ! the new velocity is ua*rdt
593
594      CALL lbc_lnk( ua, 'U', -1. )
595      CALL lbc_lnk( va, 'V', -1. )
596
597      ub(:,:,:) = ua(:,:,:) * rdt
598      vb(:,:,:) = va(:,:,:) * rdt
599      ua(:,:,:) = 0.e0
600      va(:,:,:) = 0.e0
601      un(:,:,:) = ub(:,:,:)
602      vn(:,:,:) = vb(:,:,:)
603       
604#if defined key_dynspg_rl
605      IF( lk_isl )   bsfb(:,:) = bsfn(:,:)          ! Put bsfb to zero
606#endif
607
608      ! Compute the divergence and curl
609
610      CALL div_cur( nit000 )            ! now horizontal divergence and curl
611
612      hdivb(:,:,:) = hdivn(:,:,:)       ! set the before to the now value
613      rotb (:,:,:) = rotn (:,:,:)       ! set the before to the now value
614
615   END SUBROUTINE istate_uvg
616
617   !!=====================================================================
618END MODULE istate
Note: See TracBrowser for help on using the repository browser.