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

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

nemo_v1_update_044 : CT : update the passive tracers TOP component and the standard GYRE configuration

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 28.7 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_rl
81      ! rigid-lid formulation
82      bsfb(:,:) = 0.e0      ! before barotropic stream-function
83      bsfn(:,:) = 0.e0      ! now    barotropic stream-function
84      bsfd(:,:) = 0.e0      ! barotropic stream-function trend
85#endif
86      ! free surface formulation
87      sshb(:,:) = 0.e0      ! before sea-surface height
88      sshn(:,:) = 0.e0      ! now    sea-surface height
89
90
91      IF( ln_rstart ) THEN                    ! Restart from a file
92         !                                    ! -------------------
93         neuler = 1                              ! Set time-step indicator at nit000 (leap-frog)
94         CALL rst_read                           ! Read the restart file
95      ELSE
96         !                                    ! Start from rest
97         !                                    ! ---------------
98         neuler = 0                              ! Set time-step indicator at nit000 (euler forward)
99         adatrj = 0._wp
100         IF( cp_cfg == 'eel' ) THEN
101            CALL istate_eel                      ! EEL   configuration : start from pre-defined
102            !                                    !                       velocity and thermohaline fields
103         ELSEIF( cp_cfg == 'gyre' ) THEN         
104            CALL istate_gyre                     ! GYRE  configuration : start from pre-defined temperature
105            !                                    !                       and salinity fields
106         ELSE
107         !                                       ! Other configurations: Initial temperature and salinity fields
108#if defined key_dtatem
109            CALL dta_tem( nit000 )                  ! read 3D temperature data
110            tb(:,:,:) = t_dta(:,:,:)                ! use temperature data read
111            tn(:,:,:) = t_dta(:,:,:)
112#else
113            IF(lwp) WRITE(numout,*)                 ! analytical temperature profile
114            IF(lwp) WRITE(numout,*)' Temperature initialization using an analytic profile'
115            CALL istate_tem
116#endif
117#if defined key_dtasal
118            CALL dta_sal( nit000 )                  ! read 3D salinity data
119            sb(:,:,:) = s_dta(:,:,:)                ! use salinity data read
120            sn(:,:,:) = s_dta(:,:,:)
121#else
122            ! No salinity data
123            IF(lwp)WRITE(numout,*)                  ! analytical salinity profile
124            IF(lwp)WRITE(numout,*)' Salinity initialisation using a constant value'
125            CALL istate_sal
126#endif
127         ENDIF
128
129      ENDIF
130      !                                       ! Vertical velocity
131      !                                       ! -----------------
132      CALL wzv( nit000 )                         ! from horizontal divergence
133
134   END SUBROUTINE istate_init
135
136
137   SUBROUTINE istate_tem
138      !!---------------------------------------------------------------------
139      !!                  ***  ROUTINE istate_tem  ***
140      !!   
141      !! ** Purpose :   Intialization of the temperature field with an
142      !!      analytical profile or a file (i.e. in EEL configuration)
143      !!
144      !! ** Method  :   Use Philander analytic profile of temperature
145      !!
146      !! References :  Philander ???
147      !!
148      !! History :
149      !!   4.0  !  89-12  (P. Andrich)  Original code
150      !!   6.0  !  96-01  (G. Madec)  terrain following coordinates
151      !!   9.0  !  02-09  (G. Madec)  F90: Free form
152      !!----------------------------------------------------------------------
153      !! * Local declarations
154      INTEGER :: ji, jj, jk
155      !!----------------------------------------------------------------------
156
157      IF(lwp) WRITE(numout,*)
158      IF(lwp) WRITE(numout,*) 'istate_tem : initial temperature profile'
159      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
160
161      DO jk = 1, jpk
162         DO jj = 1, jpj
163            DO ji = 1, jpi
164               tn(ji,jj,jk) = (  ( ( 7.5 - 0.*ABS(gphit(ji,jj))/30. )   &
165                  &               *( 1.-TANH((fsdept(ji,jj,jk)-80.)/30.) )   &
166                  &            + 10.*(5000.-fsdept(ji,jj,jk))/5000.)  ) * tmask(ji,jj,jk)
167               tb(ji,jj,jk) = tn(ji,jj,jk)
168          END DO
169        END DO
170      END DO
171
172      IF(lwp) CALL prizre( tn    , jpi   , jpj   , jpk   , jpj/2 ,   &
173         &                 1     , jpi   , 5     , 1     , jpk   ,   &
174         &                 1     , 1.    , numout                  )
175
176   END SUBROUTINE istate_tem
177
178
179   SUBROUTINE istate_sal
180      !!---------------------------------------------------------------------
181      !!                  ***  ROUTINE istate_sal  ***
182      !!
183      !! ** Purpose :   Intialize the salinity field with an analytic profile
184      !!
185      !! ** Method  :   Use to a constant value 35.5
186      !!             
187      !! ** Action  :   Initialize sn and sb
188      !!
189      !! History :
190      !!   4.0  !  89-12  (P. Andrich)  Original code
191      !!   8.5  !  02-09  (G. Madec)  F90: Free form
192      !!----------------------------------------------------------------------
193      !! * Local declarations
194      REAL(wp) ::   zsal = 35.50_wp
195      !!----------------------------------------------------------------------
196
197      IF(lwp) WRITE(numout,*)
198      IF(lwp) WRITE(numout,*) 'istate_sal : initial salinity : ', zsal
199      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
200
201      sn(:,:,:) = zsal * tmask(:,:,:)
202      sb(:,:,:) = sn(:,:,:)
203     
204   END SUBROUTINE istate_sal
205
206
207   SUBROUTINE istate_eel
208      !!----------------------------------------------------------------------
209      !!                   ***  ROUTINE istate_eel  ***
210      !!
211      !! ** Purpose :   Initialization of the dynamics and tracers for EEL R5
212      !!      configuration (channel with or without a topographic bump)
213      !!
214      !! ** Method  : - set temprature field
215      !!              - set salinity field
216      !!              - set velocity field including horizontal divergence
217      !!                and relative vorticity fields
218      !!
219      !! History :
220      !!   8.0  !  01-09  (M. Levy, M. Ben Jelloul)  read file for EEL 2 & 6
221      !!   9.0  !  03-09  (G. Madec, C. Talandier)  EEL 5
222      !!   9.0  !  05-11  (V. Garnier) Surface pressure gradient organization
223      !!----------------------------------------------------------------------
224      !! * Modules used
225      USE eosbn2     ! eq. of state, Brunt Vaisala frequency (eos     routine)
226      USE divcur     ! hor. divergence & rel. vorticity      (div_cur routine)
227      USE ioipsl
228 
229      !! * Local declarations
230      LOGICAL :: llog
231      CHARACTER (len=21) ::   &
232         clname = 'eel.initemp',   &  ! filename (for EEL R2 or R6)
233         clvar  = 'initemp'           ! variable name
234      INTEGER  ::   inum              ! temporary logical unit
235      INTEGER  ::   ji, jj, jk        ! dummy loop indices
236      INTEGER  ::   ilev, itime       ! temporary integers
237      REAL(wp) ::   &
238         zh1, zh2, zslope, zcst       ! temporary scalars
239      REAL(wp) ::   &
240         zt1  = 12._wp,            &  ! surface temperature value (EEL R5)
241         zt2  =  2._wp,            &  ! bottom  temperature value (EEL R5)
242         zsal = 35.5_wp                ! constant salinity (EEL R2, R5 and R6)
243      REAL(wp) ::   &
244         zdt,  zdate0                 ! temporary scalars
245      REAL(wp), DIMENSION(jpk) ::  &
246         zdept                        ! temporary workspace
247      REAL(wp), DIMENSION(jpiglo,jpjglo) ::   &
248         zlamt, zphit                 ! temporary workspace
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(  1  )
267            zh2 = gdept(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            ub(:,:,:) = 0.1 * umask(:,:,:)
297            un(:,:,:) = ub(:,:,:)
298
299            DO jj = 1, jpjglo
300               zssh(:,jj) = ( .22 - ( FLOAT(jj-3) * (0.44) ) / 99. )
301            END DO
302            DO jj = 1, nlcj
303               DO ji = 1, nlci
304                  sshb(ji,jj) = zssh( mig(ji) , mjg(jj) ) * tmask(ji,jj,1)
305               END DO
306            END DO
307            sshb(nlci+1:jpi,      :   ) = 0.e0      ! set to zero extra mpp columns
308            sshb(      :   ,nlcj+1:jpj) = 0.e0      ! set to zero extra mpp rows
309
310            sshn(:,:) = sshb(:,:)                   ! set now ssh to the before value
311
312            ! horizontal divergence and relative vorticity (curl)
313            CALL div_cur( nit000 )
314
315            ! N.B. the vertical velocity will be computed from the horizontal divergence field
316            ! in istate by a call to wzv routine
317# endif
318
319
320            !                                     ! ==========================
321         CASE ( 2 , 6 )                           ! EEL R2 or R6 configuration
322            !                                     ! ==========================
323 
324            ! set temperature field with a NetCDF file
325            ! ----------------------------------------
326            IF(lwp) WRITE(numout,*)
327            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R2 or R6: read initial temperature in a NetCDF file'
328            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
329
330            itime  = 0
331            clname = 'eel.initemp'
332#if defined key_agrif
333      if ( .NOT. Agrif_Root() ) then
334         clname = TRIM(Agrif_CFixed())//'_'//TRIM(clname)
335      endif
336#endif               
337            llog   = .FALSE.
338            ilev   = jpk
339            zlamt(:,:) = 0.e0
340            zphit(:,:) = 0.e0
341            CALL restini( clname, jpidta, jpjdta, zlamt , zphit ,   &
342               &          ilev  , zdept , clname, itime , zdate0,   &
343               &          zdt   , inum  , domain_id=nidom )
344            CALL restget( inum  , 'initemp', jpi, jpj, jpk,   &
345               &          0     , llog     , tb             )     ! read before temprature (tb)
346            CALL restclo( inum )
347 
348            tn(:,:,:) = tb(:,:,:)                                 ! set nox temperature to tb
349
350            IF(lwp) WRITE(numout,*) '               file name : ', clname
351            IF(lwp) CALL prizre( tn    , jpi   , jpj   , jpk   , jpj/2 ,   &
352               &                 1     , jpi   , 5     , 1     , jpk   ,   &
353               &                 1     , 1.    , numout                  )
354
355
356            ! set salinity field to a constant value
357            ! --------------------------------------
358            IF(lwp) WRITE(numout,*)
359            IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal
360            IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
361 
362            sn(:,:,:) = zsal * tmask(:,:,:)
363            sb(:,:,:) = sn(:,:,:)
364
365
366            IF( lk_isl ) THEN
367               ! Horizontal velocity : start from geostrophy (EEL config)
368               CALL eos( tn, sn, rhd )     ! now in situ density
369               CALL istate_uvg             ! compute geostrophic velocity
370
371               ! N.B. the vertical velocity will be computed from the horizontal divergence field
372               ! in istate by a call to wzv routine
373            ENDIF
374            !                                    ! ===========================
375         CASE DEFAULT                            ! NONE existing configuration
376            !                                    ! ===========================
377            IF(lwp) WRITE(numout,cform_err) 
378            IF(lwp) WRITE(numout,*) 'EEL with a ', jp_cfg,' km resolution is not coded'
379            nstop = nstop +1
380      END SELECT
381
382   END SUBROUTINE istate_eel
383
384
385   SUBROUTINE istate_gyre
386      !!----------------------------------------------------------------------
387      !!                   ***  ROUTINE istate_gyre  ***
388      !!
389      !! ** Purpose :   Initialization of the dynamics and tracers for GYRE
390      !!      configuration (double gyre with rotated domain)
391      !!
392      !! ** Method  : - set temprature field
393      !!              - set salinity field
394      !!
395      !! ** History :                                     
396      !!      9.0  !  04-05  (A. Koch-Larrouy)  Original code
397      !!----------------------------------------------------------------------
398      !! * Modules used
399      USE ioipsl
400
401      !! * Local variables
402      INTEGER, PARAMETER ::   jpmois = 12
403      INTEGER, PARAMETER ::   &
404         ntsinit = 0         ! (0/1) (analytical/input data files) T&S initialization
405
406      CHARACTER (len=32) ::   clname
407      INTEGER :: ji, jj, jk                       ! dummy loop indices
408      INTEGER ::   ipi, ipj, ipk, itime           ! temporary integers
409      INTEGER, DIMENSION(jpmois) ::   istep
410
411      REAL(wp) ::   zdate0, zdt
412      REAL(wp), DIMENSION(jpk) ::   zlev
413      REAL(wp), DIMENSION(jpi,jpj) ::   zlon, zlat
414      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zt_dta, zs_dta
415      !!----------------------------------------------------------------------
416
417      SELECT CASE ( ntsinit)
418
419      CASE ( 0 )                  ! analytical T/S profil deduced from LEVITUS
420         IF(lwp) WRITE(numout,*)
421         IF(lwp) WRITE(numout,*) 'istate_gyre : initial analytical T and S profil deduced from LEVITUS '
422         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
423
424         DO jk = 1, jpk
425            DO jj = 1, jpj
426               DO ji = 1, jpi
427                  tn(ji,jj,jk) = (  16. - 12. * TANH( (fsdept(ji,jj,jk) - 400) / 700 )         )   &
428                       &           * (-TANH( (500-fsdept(ji,jj,jk)) / 150 ) + 1) / 2               &
429                       &       + (      15. * ( 1. - TANH( (fsdept(ji,jj,jk)-50.) / 1500.) )       &
430                       &                - 1.4 * TANH((fsdept(ji,jj,jk)-100.) / 100.)               &   
431                       &                + 7.  * (1500. - fsdept(ji,jj,jk)) / 1500.             )   & 
432                       &           * (-TANH( (fsdept(ji,jj,jk) - 500) / 150) + 1) / 2
433                  tn(ji,jj,jk) = tn(ji,jj,jk) * tmask(ji,jj,jk)
434                  tb(ji,jj,jk) = tn(ji,jj,jk)
435
436                  sn(ji,jj,jk) =  (  36.25 - 1.13 * TANH( (fsdept(ji,jj,jk) - 305) / 460 )  )  &
437                     &              * (-TANH((500 - fsdept(ji,jj,jk)) / 150) + 1) / 2          &
438                     &          + (  35.55 + 1.25 * (5000. - fsdept(ji,jj,jk)) / 5000.         &
439                     &                - 1.62 * TANH( (fsdept(ji,jj,jk) - 60.  ) / 650. )       &
440                     &                + 0.2  * TANH( (fsdept(ji,jj,jk) - 35.  ) / 100. )       &
441                     &                + 0.2  * TANH( (fsdept(ji,jj,jk) - 1000.) / 5000.)    )  &
442                     &              * (-TANH((fsdept(ji,jj,jk) - 500) / 150) + 1) / 2 
443                  sn(ji,jj,jk) = sn(ji,jj,jk) * tmask(ji,jj,jk)
444                  sb(ji,jj,jk) = sn(ji,jj,jk)
445               END DO
446            END DO
447         END DO
448
449      CASE ( 1 )                  ! T/S data fields read in dta_tem.nc/data_sal.nc files
450         IF(lwp) WRITE(numout,*)
451         IF(lwp) WRITE(numout,*) 'istate_gyre : initial T and S read from dta_tem.nc/data_sal.nc files'
452         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
453         IF(lwp) WRITE(numout,*) '              NetCDF FORMAT'
454
455         ! Read temperature field
456         ! ----------------------
457         ! open file
458         zdt = rdt
459         clname = 'data_tem'
460         CALL flinopen(TRIM(clname), mig(1), nlci , mjg(1),  nlcj   &
461            &    , .false.     , ipi   , ipj  , ipk   , zlon        &
462            &    , zlat        , zlev  , itime, istep , zdate0      &
463            &    , zdt         , numtdt )
464
465         ! title, dimensions and tests
466         IF( ipi /= jpidta .OR. ipj /= jpjdta .OR. ipk /= jpk ) THEN
467            IF(lwp) THEN
468               WRITE(numout,*)
469               WRITE(numout,*) 'problem with dimensions'
470               WRITE(numout,*) ' ipi ',ipi,' jpidta ',jpidta
471               WRITE(numout,*) ' ipj ',ipj,' jpjdta ',jpjdta
472               WRITE(numout,*) ' ipk ',ipk,' jpk ',jpk
473            ENDIF
474            STOP 'istate_gyre'
475         ENDIF
476         IF(lwp) WRITE(numout,*) itime,istep(1),zdate0,zdt,numtdt
477
478         
479         ! Read data
480         zt_dta(:,:,:) = 0.e0
481         CALL flinget( numtdt,'votemper',jpidta,jpjdta,jpk,1,1,   &
482            &          1,mig(1),nlci,mjg(1),nlcj,zt_dta(1:nlci,1:nlcj,1:jpk))
483
484         tn(:,:,:) = zt_dta(:,:,:)*tmask(:,:,:) 
485         tb(:,:,:) = zt_dta(:,:,:)*tmask(:,:,:) 
486
487         CALL flinclo( numtdt )
488
489         IF(lwp) WRITE(numout,*)
490         IF(lwp) WRITE(numout,*) '              read temperature data ok'
491         IF(lwp) WRITE(numout,*)
492
493         ! Read salinity field
494         ! -------------------
495         ! open file
496         zdt = rdt
497         clname = 'data_sal'
498         CALL flinopen(TRIM(clname), mig(1), nlci , mjg(1),  nlcj   &
499            &    , .false.     , ipi   , ipj  , ipk   , zlon        &
500            &    , zlat        , zlev  , itime, istep , zdate0      &
501            &    , zdt         , numsdt )
502
503         ! title, dimensions and tests
504
505         IF( ipi /= jpidta .OR. ipj /= jpjdta .OR. ipk /= jpk ) THEN
506             IF(lwp) THEN
507                 WRITE(numout,*)
508                 WRITE(numout,*) 'problem with dimensions'
509                 WRITE(numout,*) ' ipi ',ipi,' jpidta ',jpidta
510                 WRITE(numout,*) ' ipj ',ipj,' jpjdta ',jpjdta
511                 WRITE(numout,*) ' ipk ',ipk,' jpk ',jpk
512             ENDIF
513             STOP 'istate_gyre'
514         ENDIF
515         IF(lwp) WRITE(numout,*) itime,istep(1),zdate0,zdt,numsdt
516
517         ! Read data
518         zs_dta(:,:,:) = 0.e0
519         CALL flinget(numsdt,'vosaline',jpidta,jpjdta,jpk,1,1,   &
520            &         1,mig(1),nlci,mjg(1),nlcj,zs_dta(1:nlci,1:nlcj,1:jpk))
521
522         sn(:,:,:)  = zs_dta(:,:,:)*tmask(:,:,:) 
523         sb(:,:,:)  = zs_dta(:,:,:)*tmask(:,:,:) 
524
525         CALL flinclo( numsdt )
526
527         IF(lwp) WRITE(numout,*)
528         IF(lwp) WRITE(numout,*) '              read salinity data ok'
529         IF(lwp) WRITE(numout,*)
530
531      END SELECT
532
533      IF(lwp) THEN
534         WRITE(numout,*)
535         WRITE(numout,*) '              Initial temperature and salinity profiles:'
536         WRITE(numout, "(9x,' level   gdept   temperature   salinity   ')" )
537         WRITE(numout, "(10x, i4, 3f10.2)" ) ( jk, gdept(jk), tn(2,2,jk), sn(2,2,jk), jk = 1, jpk )
538      ENDIF
539
540     
541   END SUBROUTINE istate_gyre
542
543
544
545   SUBROUTINE istate_uvg
546      !!----------------------------------------------------------------------
547      !!                  ***  ROUTINE istate_uvg  ***
548      !!
549      !! ** Purpose :   Compute the geostrophic velocities from (tn,sn) fields
550      !!
551      !! ** Method  :   Using the hydrostatic hypothesis the now hydrostatic
552      !!      pressure is computed by integrating the in-situ density from the
553      !!      surface to the bottom.
554      !!                 p=integral [ rau*g dz ]
555      !!
556      !! History :
557      !!   8.1  !  01-09  (M. Levy, M. Ben Jelloul)  Original code
558      !!   8.5  !  02-09  (G. Madec)  F90: Free form
559      !!   9.0  !  05-11  (V. Garnier) Surface pressure gradient organization
560      !!----------------------------------------------------------------------
561      !! * Modules used
562      USE eosbn2          ! eq. of state, Brunt Vaisala frequency (eos     routine)
563      USE dynspg          ! surface pressure gradient             (dyn_spg routine)
564      USE divcur          ! hor. divergence & rel. vorticity      (div_cur routine)
565      USE lbclnk          ! ocean lateral boundary condition (or mpp link)
566
567      !! * Local declarations
568      INTEGER ::   ji, jj, jk        ! dummy loop indices
569      INTEGER ::   indic             ! ???
570      REAL(wp) ::   &
571         zmsv, zphv, zmsu, zphu,  &  ! temporary scalars
572         zalfg
573      REAL(wp), DIMENSION (jpi,jpj,jpk) ::   &
574         zprn                        ! workspace
575      !!----------------------------------------------------------------------
576
577      IF(lwp) WRITE(numout,*) 
578      IF(lwp) WRITE(numout,*) 'istate_uvg : Start from Geostrophy'
579      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
580
581      ! Compute the now hydrostatic pressure
582      ! ------------------------------------
583
584      zalfg = 0.5 * grav * rau0
585      ! Surface value
586      zprn(:,:,1) = zalfg * fse3w(:,:,1) * ( 1 + rhd(:,:,1) )
587
588      ! Vertical integration from the surface
589      DO jk = 2, jpkm1
590         zprn(:,:,jk) = zprn(:,:,jk-1)   &
591            &         + zalfg * fse3w(:,:,jk) * ( 2. + rhd(:,:,jk) + rhd(:,:,jk-1) )
592      END DO 
593
594      ! Compute geostrophic balance
595      ! ---------------------------
596
597      DO jk = 1, jpkm1
598         DO jj = 2, jpjm1
599            DO ji = fs_2, fs_jpim1   ! vertor opt.
600               zmsv = 1. / MAX(  umask(ji-1,jj+1,jk) + umask(ji  ,jj+1,jk)   &
601                               + umask(ji-1,jj  ,jk) + umask(ji  ,jj  ,jk) , 1.  )
602               zphv = ( zprn(ji  ,jj+1,jk) - zprn(ji-1,jj+1,jk) ) * umask(ji-1,jj+1,jk) / e1u(ji-1,jj+1)   &
603                    + ( zprn(ji+1,jj+1,jk) - zprn(ji  ,jj+1,jk) ) * umask(ji  ,jj+1,jk) / e1u(ji  ,jj+1)   &
604                    + ( zprn(ji  ,jj  ,jk) - zprn(ji-1,jj  ,jk) ) * umask(ji-1,jj  ,jk) / e1u(ji-1,jj  )   &
605                    + ( zprn(ji+1,jj  ,jk) - zprn(ji  ,jj  ,jk) ) * umask(ji  ,jj  ,jk) / e1u(ji  ,jj  )
606               zphv = 1. / rau0 * zphv * zmsv * vmask(ji,jj,jk)
607
608               zmsu = 1. / MAX(  vmask(ji+1,jj  ,jk) + vmask(ji  ,jj  ,jk)   &
609                               + vmask(ji+1,jj-1,jk) + vmask(ji  ,jj-1,jk) , 1.  )
610               zphu = ( zprn(ji+1,jj+1,jk) - zprn(ji+1,jj  ,jk) ) * vmask(ji+1,jj  ,jk) / e2v(ji+1,jj  )   &
611                    + ( zprn(ji  ,jj+1,jk) - zprn(ji  ,jj  ,jk) ) * vmask(ji  ,jj  ,jk) / e2v(ji  ,jj  )   &
612                    + ( zprn(ji+1,jj  ,jk) - zprn(ji+1,jj-1,jk) ) * vmask(ji+1,jj-1,jk) / e2v(ji+1,jj-1)   &
613                    + ( zprn(ji  ,jj  ,jk) - zprn(ji  ,jj-1,jk) ) * vmask(ji  ,jj-1,jk) / e2v(ji  ,jj-1)
614               zphu = 1. / rau0 * zphu * zmsu * umask(ji,jj,jk)
615
616               ! Compute the geostrophic velocities
617               un(ji,jj,jk) = -2. * zphu / ( ff(ji,jj) + ff(ji  ,jj-1) )
618               vn(ji,jj,jk) =  2. * zphv / ( ff(ji,jj) + ff(ji-1,jj  ) )
619            END DO
620         END DO
621      END DO
622
623      IF(lwp) WRITE(numout,*) '         we force to zero bottom velocity'
624
625      ! Susbtract the bottom velocity (level jpk-1 for flat bottom case)
626      ! to have a zero bottom velocity
627
628      DO jk = 1, jpkm1
629         un(:,:,jk) = ( un(:,:,jk) - un(:,:,jpkm1) ) * umask(:,:,jk)
630         vn(:,:,jk) = ( vn(:,:,jk) - vn(:,:,jpkm1) ) * vmask(:,:,jk)
631      END DO
632
633      CALL lbc_lnk( un, 'U', -1. )
634      CALL lbc_lnk( vn, 'V', -1. )
635     
636      ub(:,:,:) = un(:,:,:)
637      vb(:,:,:) = vn(:,:,:)
638     
639      ! WARNING !!!!!
640      ! after initializing u and v, we need to calculate the initial streamfunction bsf.
641      ! Otherwise, only the trend will be computed and the model will blow up (inconsistency).
642     
643      ! to do that, we call dyn_spg with a special trick:
644      ! we fill ua and va with the velocities divided by dt,
645      ! and the streamfunction will be brought to the right
646      ! value assuming the velocities have been set up in
647      ! one time step.
648      ! we then set bsfd to zero (first guess for next step
649      ! is d(psi)/dt = 0.)
650
651      !  sets up s false trend to calculate the barotropic
652      !  streamfunction.
653
654      ua(:,:,:) = ub(:,:,:) / rdt
655      va(:,:,:) = vb(:,:,:) / rdt
656
657      ! calls dyn_spg. we assume euler time step, starting from rest.
658      indic = 0
659      CALL dyn_spg( nit000, indic )       ! surface pressure gradient
660
661      ! the new velocity is ua*rdt
662
663      CALL lbc_lnk( ua, 'U', -1. )
664      CALL lbc_lnk( va, 'V', -1. )
665
666      ub(:,:,:) = ua(:,:,:) * rdt
667      vb(:,:,:) = va(:,:,:) * rdt
668      ua(:,:,:) = 0.e0
669      va(:,:,:) = 0.e0
670      un(:,:,:) = ub(:,:,:)
671      vn(:,:,:) = vb(:,:,:)
672       
673#if defined key_dynspg_rl
674      IF( lk_isl )   bsfb(:,:) = bsfn(:,:)          ! Put bsfb to zero
675#endif
676
677      ! Compute the divergence and curl
678
679      CALL div_cur( nit000 )            ! now horizontal divergence and curl
680
681      hdivb(:,:,:) = hdivn(:,:,:)       ! set the before to the now value
682      rotb (:,:,:) = rotn (:,:,:)       ! set the before to the now value
683
684   END SUBROUTINE istate_uvg
685
686   !!=====================================================================
687END MODULE istate
Note: See TracBrowser for help on using the repository browser.