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

Last change on this file since 28 was 15, checked in by opalod, 20 years ago

CT : UPDATE001 : First major NEMO update

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