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/DEV_r1784_mid_year_merge_2010/NEMO/OPA_SRC/DYN – NEMO

source: branches/DEV_r1784_mid_year_merge_2010/NEMO/OPA_SRC/DYN/dynspg_flt.F90 @ 1976

Last change on this file since 1976 was 1976, checked in by acc, 14 years ago

ticket #684 step 6: Add in changes between the head of the DEV_r1879_mpp_rep branch and the trunk@1879.

  • Property svn:eol-style set to native
  • Property svn:executable set to *
  • Property svn:keywords set to Id
File size: 19.1 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   !!----------------------------------------------------------------------
16#if defined key_dynspg_flt   ||   defined key_esopa 
17   !!----------------------------------------------------------------------
18   !!   'key_dynspg_flt'                              filtered free surface
19   !!----------------------------------------------------------------------
20   !!   dyn_spg_flt  : update the momentum trend with the surface pressure gradient in the filtered free surface case
21   !!   flt_rst      : read/write the time-splitting restart fields in the ocean restart file
22   !!----------------------------------------------------------------------
23   USE oce             ! ocean dynamics and tracers
24   USE dom_oce         ! ocean space and time domain
25   USE zdf_oce         ! ocean vertical physics
26   USE sbc_oce         ! surface boundary condition: ocean
27   USE obc_oce         ! Lateral open boundary condition
28   USE sol_oce         ! ocean elliptic solver
29   USE phycst          ! physical constants
30   USE domvvl          ! variable volume
31   USE dynadv          ! advection
32   USE solmat          ! matrix construction for elliptic solvers
33   USE solver          ! solver initialization
34   USE solpcg          ! preconditionned conjugate gradient solver
35   USE solsor          ! Successive Over-relaxation solver
36   USE obcdyn          ! ocean open boundary condition (obc_dyn routines)
37   USE obcvol          ! ocean open boundary condition (obc_vol routines)
38   USE bdy_oce         ! Unstructured open boundaries condition
39   USE bdydyn          ! Unstructured open boundaries condition (bdy_dyn routine)
40   USE bdyvol          ! Unstructured open boundaries condition (bdy_vol routine)
41   USE cla_dynspg      ! cross land advection
42   USE in_out_manager  ! I/O manager
43   USE lib_mpp         ! distributed memory computing library
44   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
45   USE prtctl          ! Print control
46   USE agrif_opa_interp
47   USE iom
48   USE restart         ! only for lrst_oce
49   USE lib_fortran
50
51   IMPLICIT NONE
52   PRIVATE
53
54   PUBLIC   dyn_spg_flt  ! routine called by step.F90
55   PUBLIC   flt_rst      ! routine called by istate.F90
56
57   !! * Substitutions
58#  include "domzgr_substitute.h90"
59#  include "vectopt_loop_substitute.h90"
60   !!----------------------------------------------------------------------
61   !! NEMO/OPA 3.2 , LOCEAN-IPSL (2009)
62   !! $Id$
63   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
64   !!----------------------------------------------------------------------
65
66CONTAINS
67
68   SUBROUTINE dyn_spg_flt( kt, kindic )
69      !!----------------------------------------------------------------------
70      !!                  ***  routine dyn_spg_flt  ***
71      !!
72      !! ** Purpose :   Compute the now trend due to the surface pressure
73      !!      gradient in case of filtered free surface formulation  and add
74      !!      it to the general trend of momentum equation.
75      !!
76      !! ** Method  :   Filtered free surface formulation. The surface
77      !!      pressure gradient is given by:
78      !!         spgu = 1/rau0 d/dx(ps) =  1/e1u di( sshn + btda )
79      !!         spgv = 1/rau0 d/dy(ps) =  1/e2v dj( sshn + btda )
80      !!      where sshn is the free surface elevation and btda is the after
81      !!      time derivative of the free surface elevation
82      !!       -1- evaluate the surface presure trend (including the addi-
83      !!      tional force) in three steps:
84      !!        a- compute the right hand side of the elliptic equation:
85      !!            gcb = 1/(e1t e2t) [ di(e2u spgu) + dj(e1v spgv) ]
86      !!         where (spgu,spgv) are given by:
87      !!            spgu = vertical sum[ e3u (ub+ 2 rdt ua ) ]
88      !!                 - grav 2 rdt hu /e1u di[sshn + emp]
89      !!            spgv = vertical sum[ e3v (vb+ 2 rdt va) ]
90      !!                 - grav 2 rdt hv /e2v dj[sshn + emp]
91      !!         and define the first guess from previous computation :
92      !!            zbtd = btda
93      !!            btda = 2 zbtd - btdb
94      !!            btdb = zbtd
95      !!        b- compute the relative accuracy to be reached by the
96      !!         iterative solver
97      !!        c- apply the solver by a call to sol... routine
98      !!       -2- compute and add the free surface pressure gradient inclu-
99      !!      ding the additional force used to stabilize the equation.
100      !!
101      !! ** Action : - Update (ua,va) with the surf. pressure gradient trend
102      !!
103      !! References : Roullet and Madec 1999, JGR.
104      !!---------------------------------------------------------------------
105      USE oce, ONLY :   zub   => ta   ! ta used as workspace
106      USE oce, ONLY :   zvb   => sa   ! ta used as workspace
107      !!
108      INTEGER, INTENT(in   ) ::   kt       ! ocean time-step index
109      INTEGER, INTENT(  out) ::   kindic   ! solver convergence flag (<0 if not converge)
110      !!                                   
111      INTEGER  ::   ji, jj, jk           ! dummy loop indices
112      REAL(wp) ::   z2dt, z2dtg          ! temporary scalars
113      REAL(wp) ::   zgcb, zbtd   !   -          -
114      REAL(wp) ::   ztdgu, ztdgv         !   -          -
115      !!----------------------------------------------------------------------
116      !
117      IF( kt == nit000 ) THEN
118         IF(lwp) WRITE(numout,*)
119         IF(lwp) WRITE(numout,*) 'dyn_spg_flt : surface pressure gradient trend'
120         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   (free surface constant volume case)'
121       
122         ! set to zero free surface specific arrays
123         spgu(:,:) = 0.e0                     ! surface pressure gradient (i-direction)
124         spgv(:,:) = 0.e0                     ! surface pressure gradient (j-direction)
125         CALL solver_init( nit000 )           ! Elliptic solver initialisation
126
127         ! read filtered free surface arrays in restart file
128         ! when using agrif, sshn, gcx have to be read in istate
129         IF(.NOT. lk_agrif)   CALL flt_rst( nit000, 'READ' )      ! read or initialize the following fields:
130         !                                                        ! gcx, gcxb
131      ENDIF
132
133      ! Local constant initialization
134      z2dt = 2. * rdt                                             ! time step: leap-frog
135      IF( neuler == 0 .AND. kt == nit000   )   z2dt = rdt         ! time step: Euler if restart from rest
136      IF( neuler == 0 .AND. kt == nit000+1 )   CALL sol_mat( kt )
137      z2dtg  = grav * z2dt
138
139      ! Evaluate the masked next velocity (effect of the additional force not included)
140      ! --------------------------------- 
141      IF( lk_vvl ) THEN          ! variable volume  (surface pressure gradient already included in dyn_hpg)
142         !
143         IF( ln_dynadv_vec ) THEN      ! vector form : applied on velocity
144            DO jk = 1, jpkm1
145               DO jj = 2, jpjm1
146                  DO ji = fs_2, fs_jpim1   ! vector opt.
147                     ua(ji,jj,jk) = (  ub(ji,jj,jk) + z2dt * ua(ji,jj,jk)  ) * umask(ji,jj,jk)
148                     va(ji,jj,jk) = (  vb(ji,jj,jk) + z2dt * va(ji,jj,jk)  ) * vmask(ji,jj,jk)
149                  END DO
150               END DO
151            END DO
152            !
153         ELSE                          ! flux form : applied on thickness weighted velocity
154            DO jk = 1, jpkm1
155               DO jj = 2, jpjm1
156                  DO ji = fs_2, fs_jpim1   ! vector opt.
157                     ua(ji,jj,jk) = (        ub(ji,jj,jk) * fse3u_b(ji,jj,jk)      &
158                        &           + z2dt * ua(ji,jj,jk) * fse3u_n(ji,jj,jk)  )   &
159                        &         / fse3u_a(ji,jj,jk) * umask(ji,jj,jk)
160                     va(ji,jj,jk) = (        vb(ji,jj,jk) * fse3v_b(ji,jj,jk)      &
161                        &           + z2dt * va(ji,jj,jk) * fse3v_n(ji,jj,jk)  )   &
162                        &         / fse3v_a(ji,jj,jk) * vmask(ji,jj,jk)
163                 END DO
164               END DO
165            END DO
166            !
167         ENDIF
168         !
169      ELSE                       ! fixed volume  (add the surface pressure gradient + unweighted time stepping)
170         !
171         DO jj = 2, jpjm1              ! Surface pressure gradient (now)
172            DO ji = fs_2, fs_jpim1   ! vector opt.
173               spgu(ji,jj) = - grav * ( sshn(ji+1,jj) - sshn(ji,jj) ) / e1u(ji,jj)
174               spgv(ji,jj) = - grav * ( sshn(ji,jj+1) - sshn(ji,jj) ) / e2v(ji,jj)
175            END DO
176         END DO
177         DO jk = 1, jpkm1              ! unweighted time stepping
178            DO jj = 2, jpjm1
179               DO ji = fs_2, fs_jpim1   ! vector opt.
180                  ua(ji,jj,jk) = (  ub(ji,jj,jk) + z2dt * ( ua(ji,jj,jk) + spgu(ji,jj) )  ) * umask(ji,jj,jk)
181                  va(ji,jj,jk) = (  vb(ji,jj,jk) + z2dt * ( va(ji,jj,jk) + spgv(ji,jj) )  ) * vmask(ji,jj,jk)
182               END DO
183            END DO
184         END DO
185         !
186      ENDIF
187
188#if defined key_obc
189      CALL obc_dyn( kt )      ! Update velocities on each open boundary with the radiation algorithm
190      CALL obc_vol( kt )      ! Correction of the barotropic componant velocity to control the volume of the system
191#endif
192#if defined key_bdy
193      ! Update velocities on unstructured boundary using the Flow Relaxation Scheme
194      CALL bdy_dyn_frs( kt )
195
196      IF (ln_bdy_vol) THEN
197      ! Correction of the barotropic component velocity to control the volume of the system
198        CALL bdy_vol( kt )
199      ENDIF
200#endif
201#if defined key_agrif
202      CALL Agrif_dyn( kt )    ! Update velocities on each coarse/fine interfaces
203#endif
204#if defined key_orca_r2
205      IF( n_cla == 1 )   CALL dyn_spg_cla( kt )      ! Cross Land Advection (update (ua,va))
206#endif
207
208      ! compute the next vertically averaged velocity (effect of the additional force not included)
209      ! ---------------------------------------------
210      DO jj = 2, jpjm1
211         DO ji = fs_2, fs_jpim1   ! vector opt.
212            spgu(ji,jj) = 0.e0
213            spgv(ji,jj) = 0.e0
214         END DO
215      END DO
216
217      ! vertical sum
218!CDIR NOLOOPCHG
219      IF( lk_vopt_loop ) THEN          ! vector opt., forced unroll
220         DO jk = 1, jpkm1
221            DO ji = 1, jpij
222               spgu(ji,1) = spgu(ji,1) + fse3u(ji,1,jk) * ua(ji,1,jk)
223               spgv(ji,1) = spgv(ji,1) + fse3v(ji,1,jk) * va(ji,1,jk)
224            END DO
225         END DO
226      ELSE                        ! No  vector opt.
227         DO jk = 1, jpkm1
228            DO jj = 2, jpjm1
229               DO ji = 2, jpim1
230                  spgu(ji,jj) = spgu(ji,jj) + fse3u(ji,jj,jk) * ua(ji,jj,jk)
231                  spgv(ji,jj) = spgv(ji,jj) + fse3v(ji,jj,jk) * va(ji,jj,jk)
232               END DO
233            END DO
234         END DO
235      ENDIF
236
237      ! transport: multiplied by the horizontal scale factor
238      DO jj = 2, jpjm1
239         DO ji = fs_2, fs_jpim1   ! vector opt.
240            spgu(ji,jj) = spgu(ji,jj) * e2u(ji,jj)
241            spgv(ji,jj) = spgv(ji,jj) * e1v(ji,jj)
242         END DO
243      END DO
244      CALL lbc_lnk( spgu, 'U', -1. )       ! lateral boundary conditions
245      CALL lbc_lnk( spgv, 'V', -1. )
246
247      IF( lk_vvl ) CALL sol_mat( kt )      ! build the matrix at kt (vvl case only)
248
249      ! Right hand side of the elliptic equation and first guess
250      ! --------------------------------------------------------
251      DO jj = 2, jpjm1
252         DO ji = fs_2, fs_jpim1   ! vector opt.
253            ! Divergence of the after vertically averaged velocity
254            zgcb =  spgu(ji,jj) - spgu(ji-1,jj)   &
255                  + spgv(ji,jj) - spgv(ji,jj-1)
256            gcb(ji,jj) = gcdprc(ji,jj) * zgcb
257            ! First guess of the after barotropic transport divergence
258            zbtd = gcx(ji,jj)
259            gcx (ji,jj) = 2. * zbtd   - gcxb(ji,jj)
260            gcxb(ji,jj) =      zbtd
261         END DO
262      END DO
263      ! applied the lateral boundary conditions
264      IF( nn_solv == 2 .AND. MAX( jpr2di, jpr2dj ) > 0 )   CALL lbc_lnk_e( gcb, c_solver_pt, 1. )   
265
266#if defined key_agrif
267      IF( .NOT. AGRIF_ROOT() ) THEN
268         ! add contribution of gradient of after barotropic transport divergence
269         IF( nbondi == -1 .OR. nbondi == 2 )   gcb(3     ,:) =   &
270            &    gcb(3     ,:) - z2dtg * z2dt * laplacu(2     ,:) * gcdprc(3     ,:) * hu(2     ,:) * e2u(2     ,:)
271         IF( nbondi ==  1 .OR. nbondi == 2 )   gcb(nlci-2,:) =   &
272            &    gcb(nlci-2,:) + z2dtg * z2dt * laplacu(nlci-2,:) * gcdprc(nlci-2,:) * hu(nlci-2,:) * e2u(nlci-2,:)
273         IF( nbondj == -1 .OR. nbondj == 2 )   gcb(:     ,3) =   &
274            &    gcb(:,3     ) - z2dtg * z2dt * laplacv(:,2     ) * gcdprc(:,3     ) * hv(:,2     ) * e1v(:,2     )
275         IF( nbondj ==  1 .OR. nbondj == 2 )   gcb(:,nlcj-2) =   &
276            &    gcb(:,nlcj-2) + z2dtg * z2dt * laplacv(:,nlcj-2) * gcdprc(:,nlcj-2) * hv(:,nlcj-2) * e1v(:,nlcj-2)
277      ENDIF
278#endif
279
280
281      ! Relative precision (computation on one processor)
282      ! ------------------
283      rnorme =0.e0
284      rnorme = GLOB_SUM( gcb(1:jpi,1:jpj) * gcdmat(1:jpi,1:jpj) * gcb(1:jpi,1:jpj) * bmask(:,:) )
285
286      epsr = eps * eps * rnorme
287      ncut = 0
288      ! if rnorme is 0, the solution is 0, the solver is not called
289      IF( rnorme == 0.e0 ) THEN
290         gcx(:,:) = 0.e0
291         res   = 0.e0
292         niter = 0
293         ncut  = 999
294      ENDIF
295
296      ! Evaluate the next transport divergence
297      ! --------------------------------------
298      !    Iterarive solver for the elliptic equation (except IF sol.=0)
299      !    (output in gcx with boundary conditions applied)
300      kindic = 0
301      IF( ncut == 0 ) THEN
302         IF    ( nn_solv == 1 ) THEN   ;   CALL sol_pcg( kindic )      ! diagonal preconditioned conjuguate gradient
303         ELSEIF( nn_solv == 2 ) THEN   ;   CALL sol_sor( kindic )      ! successive-over-relaxation
304         ENDIF
305      ENDIF
306
307      ! Transport divergence gradient multiplied by z2dt
308      ! --------------------------------------------====
309      DO jj = 2, jpjm1
310         DO ji = fs_2, fs_jpim1   ! vector opt.
311            ! trend of Transport divergence gradient
312            ztdgu = z2dtg * (gcx(ji+1,jj  ) - gcx(ji,jj) ) / e1u(ji,jj)
313            ztdgv = z2dtg * (gcx(ji  ,jj+1) - gcx(ji,jj) ) / e2v(ji,jj)
314            ! multiplied by z2dt
315#if defined key_obc
316            ! caution : grad D = 0 along open boundaries
317            spgu(ji,jj) = z2dt * ztdgu * obcumask(ji,jj)
318            spgv(ji,jj) = z2dt * ztdgv * obcvmask(ji,jj)
319#elif defined key_bdy
320            ! caution : grad D = 0 along open boundaries
321            ! Remark: The filtering force could be reduced here in the FRS zone
322            !         by multiplying spgu/spgv by (1-alpha) ?? 
323            spgu(ji,jj) = z2dt * ztdgu * bdyumask(ji,jj)
324            spgv(ji,jj) = z2dt * ztdgv * bdyvmask(ji,jj)           
325#else
326            spgu(ji,jj) = z2dt * ztdgu
327            spgv(ji,jj) = z2dt * ztdgv
328#endif
329         END DO
330      END DO
331
332#if defined key_agrif     
333      IF( .NOT. Agrif_Root() ) THEN
334         ! caution : grad D (fine) = grad D (coarse) at coarse/fine interface
335         IF( nbondi == -1 .OR. nbondi == 2 ) spgu(2     ,:) = z2dtg * z2dt * laplacu(2     ,:) * umask(2     ,:,1)
336         IF( nbondi ==  1 .OR. nbondi == 2 ) spgu(nlci-2,:) = z2dtg * z2dt * laplacu(nlci-2,:) * umask(nlci-2,:,1)
337         IF( nbondj == -1 .OR. nbondj == 2 ) spgv(:,2     ) = z2dtg * z2dt * laplacv(:,2     ) * vmask(:     ,2,1)
338         IF( nbondj ==  1 .OR. nbondj == 2 ) spgv(:,nlcj-2) = z2dtg * z2dt * laplacv(:,nlcj-2) * vmask(:,nlcj-2,1)
339      ENDIF
340#endif     
341      ! Add the trends multiplied by z2dt to the after velocity
342      ! -------------------------------------------------------
343      !     ( c a u t i o n : (ua,va) here are the after velocity not the
344      !                       trend, the leap-frog time stepping will not
345      !                       be done in dynnxt.F90 routine)
346      DO jk = 1, jpkm1
347         DO jj = 2, jpjm1
348            DO ji = fs_2, fs_jpim1   ! vector opt.
349               ua(ji,jj,jk) = ( ua(ji,jj,jk) + spgu(ji,jj) ) * umask(ji,jj,jk)
350               va(ji,jj,jk) = ( va(ji,jj,jk) + spgv(ji,jj) ) * vmask(ji,jj,jk)
351            END DO
352         END DO
353      END DO
354
355      ! write filtered free surface arrays in restart file
356      ! --------------------------------------------------
357      IF( lrst_oce ) CALL flt_rst( kt, 'WRITE' )
358      !
359   END SUBROUTINE dyn_spg_flt
360
361
362   SUBROUTINE flt_rst( kt, cdrw )
363     !!---------------------------------------------------------------------
364     !!                   ***  ROUTINE ts_rst  ***
365     !!
366     !! ** Purpose : Read or write filtered free surface arrays in restart file
367     !!----------------------------------------------------------------------
368     INTEGER         , INTENT(in) ::   kt         ! ocean time-step
369     CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag
370     !!----------------------------------------------------------------------
371
372     IF( TRIM(cdrw) == 'READ' ) THEN
373        IF( iom_varid( numror, 'gcx', ldstop = .FALSE. ) > 0 ) THEN
374! Caution : extra-hallow
375! gcx and gcxb are defined as: DIMENSION(1-jpr2di:jpi+jpr2di,1-jpr2dj:jpj+jpr2dj)
376           CALL iom_get( numror, jpdom_autoglo, 'gcx' , gcx (1:jpi,1:jpj) )
377           CALL iom_get( numror, jpdom_autoglo, 'gcxb', gcxb(1:jpi,1:jpj) )
378           IF( neuler == 0 )   gcxb(:,:) = gcx (:,:)
379        ELSE
380           gcx (:,:) = 0.e0
381           gcxb(:,:) = 0.e0
382        ENDIF
383     ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN
384! Caution : extra-hallow
385! gcx and gcxb are defined as: DIMENSION(1-jpr2di:jpi+jpr2di,1-jpr2dj:jpj+jpr2dj)
386        CALL iom_rstput( kt, nitrst, numrow, 'gcx' , gcx (1:jpi,1:jpj) )
387        CALL iom_rstput( kt, nitrst, numrow, 'gcxb', gcxb(1:jpi,1:jpj) )
388     ENDIF
389     !
390   END SUBROUTINE flt_rst
391
392#else
393   !!----------------------------------------------------------------------
394   !!   Default case :   Empty module   No standart free surface cst volume
395   !!----------------------------------------------------------------------
396CONTAINS
397   SUBROUTINE dyn_spg_flt( kt, kindic )       ! Empty routine
398      WRITE(*,*) 'dyn_spg_flt: You should not have seen this print! error?', kt, kindic
399   END SUBROUTINE dyn_spg_flt
400   SUBROUTINE flt_rst    ( kt, cdrw )         ! Empty routine
401      INTEGER         , INTENT(in) ::   kt         ! ocean time-step
402      CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag
403      WRITE(*,*) 'flt_rst: You should not have seen this print! error?', kt, cdrw
404   END SUBROUTINE flt_rst
405#endif
406   
407   !!======================================================================
408END MODULE dynspg_flt
Note: See TracBrowser for help on using the repository browser.