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 @ 247

Last change on this file since 247 was 247, checked in by opalod, 19 years ago

CL : Add CVS Header and CeCILL licence information

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