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

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

Initial revision

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