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.
dynspg_flt.F90 in branches/UKMO/dev_r5518_GO6_package/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: branches/UKMO/dev_r5518_GO6_package/NEMOGCM/NEMO/OPA_SRC/DYN/dynspg_flt.F90 @ 11101

Last change on this file since 11101 was 11101, checked in by frrh, 5 years ago

Merge changes from Met Office GMED ticket 450 to reduce unnecessary
text output from NEMO.
This output, which is typically not switchable, is rarely of interest
in normal (non-debugging) runs and simply redunantley consumes extra
file space.
Further, the presence of this text output has been shown to
significantly degrade performance of models which are run during
Met Office HPC RAID (disk) checks.
The new code introduces switches which are configurable via the
changes made in the associated Met Office MOCI ticket 399.

File size: 20.5 KB
Line 
1MODULE dynspg_flt
2   !!======================================================================
3   !!                   ***  MODULE  dynspg_flt  ***
4   !! Ocean dynamics:  surface pressure gradient trend
5   !!======================================================================
6   !! History    OPA  !  1998-05  (G. Roullet)  free surface
7   !!                 !  1998-10  (G. Madec, M. Imbard)  release 8.2
8   !!   NEMO     O.1  !  2002-08  (G. Madec)  F90: Free form and module
9   !!             -   !  2002-11  (C. Talandier, A-M Treguier) Open boundaries
10   !!            1.0  !  2004-08  (C. Talandier) New trends organization
11   !!             -   !  2005-11  (V. Garnier) Surface pressure gradient organization
12   !!            2.0  !  2006-07  (S. Masson)  distributed restart using iom
13   !!             -   !  2006-08  (J.Chanut, A.Sellar) Calls to BDY routines.
14   !!            3.2  !  2009-03  (G. Madec, M. Leclair, R. Benshila) introduce sshwzv module
15   !!            3.7  !  2014-04  (F. Roquet, G. Madec)  add some trends diag
16   !!----------------------------------------------------------------------
17#if defined key_dynspg_flt   ||   defined key_esopa 
18   !!----------------------------------------------------------------------
19   !!   'key_dynspg_flt'                              filtered free surface
20   !!----------------------------------------------------------------------
21   !!   dyn_spg_flt  : update the momentum trend with the surface pressure gradient in the filtered free surface case
22   !!   flt_rst      : read/write the time-splitting restart fields in the ocean restart file
23   !!----------------------------------------------------------------------
24   USE oce             ! ocean dynamics and tracers
25   USE dom_oce         ! ocean space and time domain
26   USE zdf_oce         ! ocean vertical physics
27   USE sbc_oce         ! surface boundary condition: ocean
28   USE bdy_oce         ! Lateral open boundary condition
29   USE sol_oce         ! ocean elliptic solver
30   USE phycst          ! physical constants
31   USE domvvl          ! variable volume
32   USE dynadv          ! advection
33   USE solmat          ! matrix construction for elliptic solvers
34   USE solpcg          ! preconditionned conjugate gradient solver
35   USE solsor          ! Successive Over-relaxation solver
36   USE bdydyn          ! ocean open boundary condition on dynamics
37   USE bdyvol          ! ocean open boundary condition (bdy_vol routine)
38   USE cla             ! cross land advection
39   USE trd_oce         ! trends: ocean variables
40   USE trddyn          ! trend manager: dynamics
41   !
42   USE in_out_manager  ! I/O manager
43   USE lib_mpp         ! distributed memory computing library
44   USE wrk_nemo        ! Memory Allocation
45   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
46   USE prtctl          ! Print control
47   USE iom
48   USE lib_fortran
49   USE timing          ! Timing
50#if defined key_agrif
51   USE agrif_opa_interp
52#endif
53
54   IMPLICIT NONE
55   PRIVATE
56
57   PUBLIC   dyn_spg_flt  ! routine called by step.F90
58   PUBLIC   flt_rst      ! routine called by istate.F90
59
60   !! * Substitutions
61#  include "domzgr_substitute.h90"
62#  include "vectopt_loop_substitute.h90"
63   !!----------------------------------------------------------------------
64   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
65   !! $Id$
66   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
67   !!----------------------------------------------------------------------
68CONTAINS
69
70   SUBROUTINE dyn_spg_flt( kt, kindic )
71      !!----------------------------------------------------------------------
72      !!                  ***  routine dyn_spg_flt  ***
73      !!
74      !! ** Purpose :   Compute the now trend due to the surface pressure
75      !!      gradient in case of filtered free surface formulation  and add
76      !!      it to the general trend of momentum equation.
77      !!
78      !! ** Method  :   Filtered free surface formulation. The surface
79      !!      pressure gradient is given by:
80      !!         spgu = 1/rau0 d/dx(ps) =  1/e1u di( sshn + btda )
81      !!         spgv = 1/rau0 d/dy(ps) =  1/e2v dj( sshn + btda )
82      !!      where sshn is the free surface elevation and btda is the after
83      !!      time derivative of the free surface elevation
84      !!       -1- evaluate the surface presure trend (including the addi-
85      !!      tional force) in three steps:
86      !!        a- compute the right hand side of the elliptic equation:
87      !!            gcb = 1/(e1t e2t) [ di(e2u spgu) + dj(e1v spgv) ]
88      !!         where (spgu,spgv) are given by:
89      !!            spgu = vertical sum[ e3u (ub+ 2 rdt ua ) ]
90      !!                 - grav 2 rdt hu /e1u di[sshn + (emp-rnf)]
91      !!            spgv = vertical sum[ e3v (vb+ 2 rdt va) ]
92      !!                 - grav 2 rdt hv /e2v dj[sshn + (emp-rnf)]
93      !!         and define the first guess from previous computation :
94      !!            zbtd = btda
95      !!            btda = 2 zbtd - btdb
96      !!            btdb = zbtd
97      !!        b- compute the relative accuracy to be reached by the
98      !!         iterative solver
99      !!        c- apply the solver by a call to sol... routine
100      !!       -2- compute and add the free surface pressure gradient inclu-
101      !!      ding the additional force used to stabilize the equation.
102      !!
103      !! ** Action : - Update (ua,va) with the surf. pressure gradient trend
104      !!
105      !! References : Roullet and Madec, JGR, 2000.
106      !!---------------------------------------------------------------------
107      INTEGER, INTENT(in   ) ::   kt       ! ocean time-step index
108      INTEGER, INTENT(  out) ::   kindic   ! solver convergence flag (<0 if not converge)
109      !
110      INTEGER  ::   ji, jj, jk   ! dummy loop indices
111      REAL(wp) ::   z2dt, z2dtg, zgcb, zbtd, ztdgu, ztdgv   ! local scalars
112      REAL(wp), POINTER, DIMENSION(:,:,:) ::  ztrdu, ztrdv
113      REAL(wp), POINTER, DIMENSION(:,:)   ::  zpw
114      !!----------------------------------------------------------------------
115      !
116      IF( nn_timing == 1 )  CALL timing_start('dyn_spg_flt')
117      !
118      IF( kt == nit000 ) THEN
119         IF(lwp) WRITE(numout,*)
120         IF(lwp) WRITE(numout,*) 'dyn_spg_flt : surface pressure gradient trend'
121         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   (free surface constant volume case)'
122         IF(lwp .AND. lflush) CALL flush(numout)
123       
124         ! set to zero free surface specific arrays
125         spgu(:,:) = 0._wp                     ! surface pressure gradient (i-direction)
126         spgv(:,:) = 0._wp                     ! surface pressure gradient (j-direction)
127
128         ! read filtered free surface arrays in restart file
129         ! when using agrif, sshn, gcx have to be read in istate
130         IF(.NOT. lk_agrif)   CALL flt_rst( nit000, 'READ' )      ! read or initialize the following fields:
131         !                                                        ! gcx, gcxb
132      ENDIF
133
134      ! Local constant initialization
135      z2dt = 2. * rdt                                             ! time step: leap-frog
136      IF( neuler == 0 .AND. kt == nit000   )   z2dt = rdt         ! time step: Euler if restart from rest
137      IF( neuler == 0 .AND. kt == nit000+1 )   CALL sol_mat( kt )
138      z2dtg  = grav * z2dt
139
140      ! Evaluate the masked next velocity (effect of the additional force not included)
141      ! --------------------------------- 
142      IF( lk_vvl ) THEN          ! variable volume  (surface pressure gradient already included in dyn_hpg)
143         !
144         IF( ln_dynadv_vec ) THEN      ! vector form : applied on velocity
145            DO jk = 1, jpkm1
146               DO jj = 2, jpjm1
147                  DO ji = fs_2, fs_jpim1   ! vector opt.
148                     ua(ji,jj,jk) = (  ub(ji,jj,jk) + z2dt * ua(ji,jj,jk)  ) * umask(ji,jj,jk)
149                     va(ji,jj,jk) = (  vb(ji,jj,jk) + z2dt * va(ji,jj,jk)  ) * vmask(ji,jj,jk)
150                  END DO
151               END DO
152            END DO
153            !
154         ELSE                          ! flux form : applied on thickness weighted velocity
155            DO jk = 1, jpkm1
156               DO jj = 2, jpjm1
157                  DO ji = fs_2, fs_jpim1   ! vector opt.
158                     ua(ji,jj,jk) = (        ub(ji,jj,jk) * fse3u_b(ji,jj,jk)      &
159                        &           + z2dt * ua(ji,jj,jk) * fse3u_n(ji,jj,jk)  )   &
160                        &         / fse3u_a(ji,jj,jk) * umask(ji,jj,jk)
161                     va(ji,jj,jk) = (        vb(ji,jj,jk) * fse3v_b(ji,jj,jk)      &
162                        &           + z2dt * va(ji,jj,jk) * fse3v_n(ji,jj,jk)  )   &
163                        &         / fse3v_a(ji,jj,jk) * vmask(ji,jj,jk)
164                 END DO
165               END DO
166            END DO
167            !
168         ENDIF
169        IF( l_trddyn )   THEN                      ! Put here so code doesn't crash when doing KE trend but needs to be done properly
170            CALL wrk_alloc( jpi, jpj, jpk, ztrdu, ztrdv )
171         ENDIF
172         !
173      ELSE                       ! fixed volume  (add the surface pressure gradient + unweighted time stepping)
174         !
175         DO jj = 2, jpjm1              ! Surface pressure gradient (now)
176            DO ji = fs_2, fs_jpim1   ! vector opt.
177               spgu(ji,jj) = - grav * ( sshn(ji+1,jj) - sshn(ji,jj) ) / e1u(ji,jj)
178               spgv(ji,jj) = - grav * ( sshn(ji,jj+1) - sshn(ji,jj) ) / e2v(ji,jj)
179            END DO
180         END DO
181         DO jk = 1, jpkm1              ! unweighted time stepping
182            DO jj = 2, jpjm1
183               DO ji = fs_2, fs_jpim1   ! vector opt.
184                  ua(ji,jj,jk) = (  ub(ji,jj,jk) + z2dt * ( ua(ji,jj,jk) + spgu(ji,jj) )  ) * umask(ji,jj,jk)
185                  va(ji,jj,jk) = (  vb(ji,jj,jk) + z2dt * ( va(ji,jj,jk) + spgv(ji,jj) )  ) * vmask(ji,jj,jk)
186               END DO
187            END DO
188         END DO
189         !
190         IF( l_trddyn )   THEN                      ! temporary save of spg trends
191            CALL wrk_alloc( jpi, jpj, jpk, ztrdu, ztrdv )
192            DO jk = 1, jpkm1              ! unweighted time stepping
193               DO jj = 2, jpjm1
194                  DO ji = fs_2, fs_jpim1   ! vector opt.
195                     ztrdu(ji,jj,jk) = spgu(ji,jj) * umask(ji,jj,jk)
196                     ztrdv(ji,jj,jk) = spgv(ji,jj) * vmask(ji,jj,jk)
197                  END DO
198               END DO
199            END DO
200            CALL trd_dyn( ztrdu, ztrdv, jpdyn_spgexp, kt )
201         ENDIF
202         !
203      ENDIF
204
205#if defined key_bdy
206      IF( lk_bdy ) CALL bdy_dyn( kt )   ! Update velocities on each open boundary
207      IF( lk_bdy ) CALL bdy_vol( kt )   ! Correction of the barotropic component velocity to control the volume of the system
208#endif
209#if defined key_agrif
210      CALL Agrif_dyn( kt )    ! Update velocities on each coarse/fine interfaces
211#endif
212      IF( nn_cla == 1 .AND. cp_cfg == 'orca' .AND. jp_cfg == 2 )   CALL cla_dynspg( kt )      ! Cross Land Advection (update (ua,va))
213
214      ! compute the next vertically averaged velocity (effect of the additional force not included)
215      ! ---------------------------------------------
216      DO jj = 2, jpjm1
217         DO ji = fs_2, fs_jpim1   ! vector opt.
218            spgu(ji,jj) = fse3u_a(ji,jj,1) * ua(ji,jj,1)
219            spgv(ji,jj) = fse3v_a(ji,jj,1) * va(ji,jj,1)
220         END DO
221      END DO
222      DO jk = 2, jpkm1                     ! vertical sum
223         DO jj = 2, jpjm1
224            DO ji = fs_2, fs_jpim1   ! vector opt.
225               spgu(ji,jj) = spgu(ji,jj) + fse3u_a(ji,jj,jk) * ua(ji,jj,jk)
226               spgv(ji,jj) = spgv(ji,jj) + fse3v_a(ji,jj,jk) * va(ji,jj,jk)
227            END DO
228         END DO
229      END DO
230
231      DO jj = 2, jpjm1                     ! transport: multiplied by the horizontal scale factor
232         DO ji = fs_2, fs_jpim1   ! vector opt.
233            spgu(ji,jj) = spgu(ji,jj) * e2u(ji,jj)
234            spgv(ji,jj) = spgv(ji,jj) * e1v(ji,jj)
235         END DO
236      END DO
237      CALL lbc_lnk( spgu, 'U', -1. )       ! lateral boundary conditions
238      CALL lbc_lnk( spgv, 'V', -1. )
239
240      IF( lk_vvl ) CALL sol_mat( kt )      ! build the matrix at kt (vvl case only)
241
242      ! Right hand side of the elliptic equation and first guess
243      ! --------------------------------------------------------
244      DO jj = 2, jpjm1
245         DO ji = fs_2, fs_jpim1   ! vector opt.
246            ! Divergence of the after vertically averaged velocity
247            zgcb =  spgu(ji,jj) - spgu(ji-1,jj)   &
248                  + spgv(ji,jj) - spgv(ji,jj-1)
249            gcb(ji,jj) = gcdprc(ji,jj) * zgcb
250            ! First guess of the after barotropic transport divergence
251            zbtd = gcx(ji,jj)
252            gcx (ji,jj) = 2. * zbtd   - gcxb(ji,jj)
253            gcxb(ji,jj) =      zbtd
254         END DO
255      END DO
256      ! applied the lateral boundary conditions
257      IF( nn_solv == 2 .AND. MAX( jpr2di, jpr2dj ) > 0 )   CALL lbc_lnk_e( gcb, c_solver_pt, 1., jpr2di, jpr2dj )   
258
259#if defined key_agrif
260      IF( .NOT. AGRIF_ROOT() ) THEN
261         ! add contribution of gradient of after barotropic transport divergence
262         IF( nbondi == -1 .OR. nbondi == 2 )   gcb(3     ,:) =   &
263            &    gcb(3     ,:) - z2dtg * z2dt * laplacu(2     ,:) * gcdprc(3     ,:) * hu(2     ,:) * e2u(2     ,:)
264         IF( nbondi ==  1 .OR. nbondi == 2 )   gcb(nlci-2,:) =   &
265            &    gcb(nlci-2,:) + z2dtg * z2dt * laplacu(nlci-2,:) * gcdprc(nlci-2,:) * hu(nlci-2,:) * e2u(nlci-2,:)
266         IF( nbondj == -1 .OR. nbondj == 2 )   gcb(:     ,3) =   &
267            &    gcb(:,3     ) - z2dtg * z2dt * laplacv(:,2     ) * gcdprc(:,3     ) * hv(:,2     ) * e1v(:,2     )
268         IF( nbondj ==  1 .OR. nbondj == 2 )   gcb(:,nlcj-2) =   &
269            &    gcb(:,nlcj-2) + z2dtg * z2dt * laplacv(:,nlcj-2) * gcdprc(:,nlcj-2) * hv(:,nlcj-2) * e1v(:,nlcj-2)
270      ENDIF
271#endif
272
273
274      ! Relative precision (computation on one processor)
275      ! ------------------
276      rnorme =0.e0
277      rnorme = GLOB_SUM( gcb(1:jpi,1:jpj) * gcdmat(1:jpi,1:jpj) * gcb(1:jpi,1:jpj) * bmask(:,:) )
278
279      epsr = eps * eps * rnorme
280      ncut = 0
281      ! if rnorme is 0, the solution is 0, the solver is not called
282      IF( rnorme == 0._wp ) THEN
283         gcx(:,:) = 0._wp
284         res   = 0._wp
285         niter = 0
286         ncut  = 999
287      ENDIF
288
289      ! Evaluate the next transport divergence
290      ! --------------------------------------
291      !    Iterarive solver for the elliptic equation (except IF sol.=0)
292      !    (output in gcx with boundary conditions applied)
293      kindic = 0
294      IF( ncut == 0 ) THEN
295         IF    ( nn_solv == 1 ) THEN   ;   CALL sol_pcg( kindic )      ! diagonal preconditioned conjuguate gradient
296         ELSEIF( nn_solv == 2 ) THEN   ;   CALL sol_sor( kindic )      ! successive-over-relaxation
297         ENDIF
298      ENDIF
299
300      ! Transport divergence gradient multiplied by z2dt
301      ! --------------------------------------------====
302      DO jj = 2, jpjm1
303         DO ji = fs_2, fs_jpim1   ! vector opt.
304            ! trend of Transport divergence gradient
305            ztdgu = z2dtg * (gcx(ji+1,jj  ) - gcx(ji,jj) ) / e1u(ji,jj)
306            ztdgv = z2dtg * (gcx(ji  ,jj+1) - gcx(ji,jj) ) / e2v(ji,jj)
307            ! multiplied by z2dt
308#if defined key_bdy
309            IF(lk_bdy) THEN
310            ! caution : grad D = 0 along open boundaries
311               spgu(ji,jj) = z2dt * ztdgu * bdyumask(ji,jj)
312               spgv(ji,jj) = z2dt * ztdgv * bdyvmask(ji,jj)
313            ELSE
314               spgu(ji,jj) = z2dt * ztdgu
315               spgv(ji,jj) = z2dt * ztdgv
316            ENDIF
317#else
318            spgu(ji,jj) = z2dt * ztdgu
319            spgv(ji,jj) = z2dt * ztdgv
320#endif
321         END DO
322      END DO
323
324#if defined key_agrif     
325      IF( .NOT. Agrif_Root() ) THEN
326         ! caution : grad D (fine) = grad D (coarse) at coarse/fine interface
327         IF( nbondi == -1 .OR. nbondi == 2 ) spgu(2     ,:) = z2dtg * z2dt * laplacu(2     ,:) * umask(2     ,:,1)
328         IF( nbondi ==  1 .OR. nbondi == 2 ) spgu(nlci-2,:) = z2dtg * z2dt * laplacu(nlci-2,:) * umask(nlci-2,:,1)
329         IF( nbondj == -1 .OR. nbondj == 2 ) spgv(:,2     ) = z2dtg * z2dt * laplacv(:,2     ) * vmask(:     ,2,1)
330         IF( nbondj ==  1 .OR. nbondj == 2 ) spgv(:,nlcj-2) = z2dtg * z2dt * laplacv(:,nlcj-2) * vmask(:,nlcj-2,1)
331      ENDIF
332#endif     
333
334      IF( l_trddyn )   THEN                     
335         ztrdu(:,:,:) = ua(:,:,:)                 ! save the after velocity before the filtered SPG
336         ztrdv(:,:,:) = va(:,:,:)
337         !
338         CALL wrk_alloc( jpi, jpj, zpw )
339         !
340         zpw(:,:) = - z2dt * gcx(:,:)
341         CALL iom_put( "ssh_flt" , zpw )          ! output equivalent ssh modification due to implicit filter
342         !
343         !                                        ! save surface pressure flux: -pw at z=0
344         zpw(:,:) = - rau0 * grav * sshn(:,:) * wn(:,:,1) * tmask(:,:,1)
345         CALL iom_put( "pw0_exp" , zpw )
346         zpw(:,:) = wn(:,:,1)
347         CALL iom_put( "w0" , zpw )
348         zpw(:,:) =  rau0 * z2dtg * gcx(:,:) * wn(:,:,1) * tmask(:,:,1)
349         CALL iom_put( "pw0_flt" , zpw )
350         !
351         CALL wrk_dealloc( jpi, jpj, zpw ) 
352         !                                   
353      ENDIF
354     
355      ! Add the trends multiplied by z2dt to the after velocity
356      ! -------------------------------------------------------
357      !     ( c a u t i o n : (ua,va) here are the after velocity not the
358      !                       trend, the leap-frog time stepping will not
359      !                       be done in dynnxt.F90 routine)
360      DO jk = 1, jpkm1
361         DO jj = 2, jpjm1
362            DO ji = fs_2, fs_jpim1   ! vector opt.
363               ua(ji,jj,jk) = ( ua(ji,jj,jk) + spgu(ji,jj) ) * umask(ji,jj,jk)
364               va(ji,jj,jk) = ( va(ji,jj,jk) + spgv(ji,jj) ) * vmask(ji,jj,jk)
365            END DO
366         END DO
367      END DO
368
369      IF( l_trddyn )   THEN                      ! save the explicit SPG trends for further diagnostics
370         ztrdu(:,:,:) = ( ua(:,:,:) - ztrdu(:,:,:) ) / z2dt
371         ztrdv(:,:,:) = ( va(:,:,:) - ztrdv(:,:,:) ) / z2dt
372         CALL trd_dyn( ztrdu, ztrdv, jpdyn_spgflt, kt )
373         !
374         CALL wrk_dealloc( jpi, jpj, jpk, ztrdu, ztrdv ) 
375      ENDIF
376
377      IF( lrst_oce )   CALL flt_rst( kt, 'WRITE' )      ! write filtered free surface arrays in restart file
378      !
379      IF( nn_timing == 1 )   CALL timing_stop('dyn_spg_flt')
380      !
381   END SUBROUTINE dyn_spg_flt
382
383
384   SUBROUTINE flt_rst( kt, cdrw )
385      !!---------------------------------------------------------------------
386      !!                   ***  ROUTINE ts_rst  ***
387      !!
388      !! ** Purpose : Read or write filtered free surface arrays in restart file
389      !!----------------------------------------------------------------------
390      INTEGER         , INTENT(in) ::   kt     ! ocean time-step
391      CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
392      !!----------------------------------------------------------------------
393      !
394      IF( TRIM(cdrw) == 'READ' ) THEN
395         IF( iom_varid( numror, 'gcx', ldstop = .FALSE. ) > 0 ) THEN
396! Caution : extra-hallow
397! gcx and gcxb are defined as: DIMENSION(1-jpr2di:jpi+jpr2di,1-jpr2dj:jpj+jpr2dj)
398            IF(nn_timing == 2)  CALL timing_start('iom_rstget')
399            CALL iom_get( numror, jpdom_autoglo, 'gcx' , gcx (1:jpi,1:jpj) )
400            CALL iom_get( numror, jpdom_autoglo, 'gcxb', gcxb(1:jpi,1:jpj) )
401            IF(nn_timing == 2)  CALL timing_stop('iom_rstget')
402            IF( neuler == 0 )   gcxb(:,:) = gcx (:,:)
403         ELSE
404            gcx (:,:) = 0.e0
405            gcxb(:,:) = 0.e0
406         ENDIF
407      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN
408! Caution : extra-hallow
409! gcx and gcxb are defined as: DIMENSION(1-jpr2di:jpi+jpr2di,1-jpr2dj:jpj+jpr2dj)
410         IF(nn_timing == 2)  CALL timing_start('iom_rstput')
411         CALL iom_rstput( kt, nitrst, numrow, 'gcx' , gcx (1:jpi,1:jpj) )
412         CALL iom_rstput( kt, nitrst, numrow, 'gcxb', gcxb(1:jpi,1:jpj) )
413         IF(nn_timing == 2)  CALL timing_stop('iom_rstput')
414      ENDIF
415      !
416   END SUBROUTINE flt_rst
417
418#else
419   !!----------------------------------------------------------------------
420   !!   Default case :   Empty module   No standart free surface cst volume
421   !!----------------------------------------------------------------------
422CONTAINS
423   SUBROUTINE dyn_spg_flt( kt, kindic )       ! Empty routine
424      WRITE(*,*) 'dyn_spg_flt: You should not have seen this print! error?', kt, kindic
425   END SUBROUTINE dyn_spg_flt
426   SUBROUTINE flt_rst    ( kt, cdrw )         ! Empty routine
427      INTEGER         , INTENT(in) ::   kt         ! ocean time-step
428      CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag
429      WRITE(*,*) 'flt_rst: You should not have seen this print! error?', kt, cdrw
430   END SUBROUTINE flt_rst
431#endif
432   
433   !!======================================================================
434END MODULE dynspg_flt
Note: See TracBrowser for help on using the repository browser.