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_fsc_atsk.F90 in trunk/NEMO/OPA_SRC/DYN – NEMO

source: trunk/NEMO/OPA_SRC/DYN/dynspg_fsc_atsk.F90 @ 247

Last change on this file since 247 was 247, checked in by opalod, 19 years ago

CL : Add CVS Header and CeCILL licence information

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 18.2 KB
Line 
1MODULE dynspg_fsc_atsk
2   !!======================================================================
3   !!                  ***  MODULE  dynspg_fsc_atsk  ***
4   !! Ocean dynamics:  surface pressure gradient trend
5   !!======================================================================
6#if ( defined key_dynspg_fsc && defined key_autotasking )   ||   defined key_esopa
7   !!----------------------------------------------------------------------
8   !!   'key_dynspg_fsc' & 'key_autotasking'   free surface cst volume
9   !!                                                   j-k-i loop (j-slab)
10   !!----------------------------------------------------------------------
11   !!   dyn_spg_fsc_atsk : Update the momentum trend with the surface pressure
12   !!                      gradient for the free surf. constant volume case
13   !!                      with auto-tasking (j-slab) (no vectior opt.)
14   !!----------------------------------------------------------------------
15   !!   OPA 9.0 , LOCEAN-IPSL (2005)
16   !! $Header$
17   !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt
18   !!----------------------------------------------------------------------
19   !! * Modules used
20   USE oce             ! ocean dynamics and tracers
21   USE dom_oce         ! ocean space and time domain
22   USE zdf_oce         ! ocean vertical physics
23   USE trdmod          ! ocean dynamics trends
24   USE trdmod_oce      ! ocean variables trends
25   USE in_out_manager  ! I/O manager
26   USE phycst          ! physical constant
27   USE ocesbc          ! Ocean Surface Boundary condition
28   USE flxrnf          ! ocean runoffs
29   USE sol_oce         ! ocean elliptic solver
30   USE solpcg          ! preconditionned conjugate gradient solver
31   USE solsor          ! Successive Over-relaxation solver
32   USE solfet          ! FETI solver
33   USE obc_oce         ! Lateral open boundary condition
34   USE obcdyn          ! ocean open boundary condition (obc_dyn routines)
35   USE obcvol          ! ocean open boundary condition (obc_vol routines)
36   USE lib_mpp         ! distributed memory computing library
37   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
38   USE cla_dynspg      ! cross land advection
39
40   IMPLICIT NONE
41   PRIVATE
42
43   !! * Accessibility
44   PUBLIC dyn_spg_fsc_atsk        ! routine called by step.F90
45
46   !! * Shares module variables
47   LOGICAL, PUBLIC, PARAMETER ::   lk_dynspg_fsc_tsk = .TRUE.    !: free surf. cst vol. flag
48
49   !! * Substitutions
50#  include "domzgr_substitute.h90"
51   !!----------------------------------------------------------------------
52   !!   OPA 9.0 , LOCEAN-IPSL (2005)
53   !! $Header$
54   !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt
55   !!----------------------------------------------------------------------
56
57CONTAINS
58
59   SUBROUTINE dyn_spg_fsc_atsk( kt, kindic )
60      !!----------------------------------------------------------------------
61      !!                  ***  routine dyn_spg_fsc_atsk  ***
62      !!
63      !! ** Purpose :   Compute the now trend due to the surface pressure
64      !!      gradient for free surface formulation with a constant ocean
65      !!      volume case, add it to the general trend of momentum equation.
66      !!
67      !! ** Method  :   Free surface formulation. The surface pressure gradient
68      !!      is given by:
69      !!         spgu = 1/rau0 d/dx(ps) =  1/e1u di( etn + rnu btda )
70      !!         spgv = 1/rau0 d/dy(ps) =  1/e2v dj( etn + rnu btda )
71      !!      where etn is the free surface elevation and btda is the after
72      !!      of the free surface elevation
73      !!       -1- compute the after sea surface elevation from the cinematic
74      !!      surface boundary condition:
75      !!              zssha = sshb + 2 rdt ( wn - emp )
76      !!           Time filter applied on now sea surface elevation to avoid
77      !!      the divergence of two consecutive time-steps and swap of free
78      !!      surface arrays to start the next time step:
79      !!              sshb = sshn + atfp * [ sshb + zssha - 2 sshn ]
80      !!              sshn = zssha
81      !!       -2- evaluate the surface presure trend (including the addi-
82      !!      tional force) in three steps:
83      !!        a- compute the right hand side of the elliptic equation:
84      !!            gcb = 1/(e1t e2t) [ di(e2u spgu) + dj(e1v spgv) ]
85      !!         where (spgu,spgv) are given by:
86      !!            spgu = vertical sum[ e3u (ub+ 2 rdt ua ) ]
87      !!                 - grav 2 rdt hu /e1u di[sshn + emp]
88      !!            spgv = vertical sum[ e3v (vb+ 2 rdt va) ]
89      !!                 - grav 2 rdt hv /e2v dj[sshn + emp]
90      !!         and define the first guess from previous computation :
91      !!            zbtd = btda
92      !!            btda = 2 zbtd - btdb
93      !!            btdb = zbtd
94      !!        b- compute the relative accuracy to be reached by the
95      !!         iterative solver
96      !!        c- apply the solver by a call to sol... routine
97      !!       -3- compute and add the free surface pressure gradient inclu-
98      !!      ding the additional force used to stabilize the equation.
99      !!      several slabs used ('key-autotasking')
100      !!
101      !! ** Action : - Update (ua,va) with the surf. pressure gradient trend
102      !!             - Save the trends in (ztdua,ztdva) ('key_trddyn')
103      !!
104      !! References :
105      !!      Roullet and Madec 1999, JGR.
106      !!
107      !! History :
108      !!        !  98-05 (G. Roullet)  Original code
109      !!        !  98-10 (G. Madec, M. Imbard)  release 8.2
110      !!   8.5  !  02-08 (G. Madec)  F90: Free form and module
111      !!        !  02-11 (C. Talandier, A-M Treguier) Open boundaries
112      !!   9.0  !  04-08 (C. Talandier) New trends organization
113      !!---------------------------------------------------------------------
114      !! * Modules used     
115      USE oce, ONLY :    ztdua => ta,      & ! use ta as 3D workspace   
116                         ztdva => sa         ! use sa as 3D workspace   
117
118      !! * Arguments
119      INTEGER, INTENT( in ) ::   kt          ! ocean time-step index
120      INTEGER, INTENT( out ) ::   kindic     ! solver convergence flag
121                                             ! if the solver doesn t converge
122                                             ! the flag is < 0
123      !! * Local declarations
124      INTEGER  ::   ji, jj, jk               ! dummy loop indices
125      REAL(wp) ::   &              ! temporary scalars
126         z2dt, z2dtg, zraur, znugdt, znurau,   &
127         zssha, zspgu, zspgv,                  &
128         zgcb, zbtd, ztdgu, ztdgv, zgwgt
129      !!----------------------------------------------------------------------
130
131      IF( kt == nit000 ) THEN
132         IF(lwp) WRITE(numout,*)
133         IF(lwp) WRITE(numout,*) 'dyn_spg_fsc_atsk : surface pressure gradient trend'
134         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~~~~  (free surface constant volume, autotasking case)'
135
136         ! set to zero free surface specific arrays
137         spgu(:,:) = 0.e0      ! surface pressure gradient (i-direction)
138         spgv(:,:) = 0.e0      ! surface pressure gradient (j-direction)
139      ENDIF
140
141      ! 0. Local constant initialization
142      ! --------------------------------
143      ! time step: leap-frog
144      z2dt = 2. * rdt
145      ! time step: Euler if restart from rest
146      IF( neuler == 0 .AND. kt == nit000 ) z2dt = rdt
147      ! coefficients
148      z2dtg  = grav * z2dt
149      zraur  = 1. / rauw
150      znugdt =  rnu * grav * z2dt
151      znurau =  znugdt * zraur
152      IF( lk_mpp ) THEN
153         ! Mpp : export boundary values of to neighboring processors
154         !!bug ???  why only in mpp?  is it really needed???
155         CALL lbc_lnk( ua, 'U' , -1. )
156         CALL lbc_lnk( va, 'V' , -1. )
157      ENDIF
158
159      !                                                ! ===============
160      DO jj = 2, jpjm1                                 !  Vertical slab
161         !                                             ! ===============
162         ! 1. Surface pressure gradient (now)
163         ! ----------------------------
164         DO ji = 2, jpim1
165            zspgu =   - grav * ( sshn(ji+1,jj) - sshn(ji,jj) ) / e1u(ji,jj)
166            zspgv =   - grav * ( sshn(ji,jj+1) - sshn(ji,jj) ) / e2v(ji,jj)
167            spgu(ji,jj) = zspgu
168            spgv(ji,jj) = zspgv
169         END DO 
170
171         ! 2. Add the surface pressure trend to the general trend
172         ! ------------------------------------------------------
173         DO jk = 1, jpkm1
174            DO ji = 2, jpim1
175               ua(ji,jj,jk) = ua(ji,jj,jk) + spgu(ji,jj)
176               va(ji,jj,jk) = va(ji,jj,jk) + spgv(ji,jj)
177            END DO
178         END DO
179         !                                             ! ===============
180      END DO                                           !   End of slab
181      !                                                ! ===============
182
183      ! Save the surface pressure gradient trend for diagnostics
184      IF( l_trddyn )   THEN
185         DO jk = 1, jpkm1
186            ztdua(:,:,jk) = spgu(:,:) 
187            ztdva(:,:,jk) = spgv(:,:) 
188         END DO
189      ENDIF
190
191      !,,,,,,,,,,,,,,,,,,,,,,,,,,,,,synchro,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,
192
193      !                                                ! ===============
194      DO jj = 2, jpjm1                                 !  Vertical slab
195         !                                             ! ===============
196
197         ! 1. Evaluate the masked next velocity
198         ! ------------------------------------
199         !     (effect of the additional force not included)
200
201         DO jk = 1, jpkm1
202            DO ji = 2, jpim1
203               ua(ji,jj,jk) = ( ub(ji,jj,jk) + z2dt * ua(ji,jj,jk) ) * umask(ji,jj,jk)
204               va(ji,jj,jk) = ( vb(ji,jj,jk) + z2dt * va(ji,jj,jk) ) * vmask(ji,jj,jk)
205            END DO
206         END DO
207
208         !                                             ! ===============
209      END DO                                           !   End of slab
210      !                                                ! ===============
211
212#if defined key_obc
213         ! Update velocities on each open boundary with the radiation algorithm
214         CALL obc_dyn( kt )
215         ! Correction of the barotropic componant velocity to control the volume of the system
216         CALL obc_vol( kt )
217#endif
218#if defined key_orca_r2
219      IF( n_cla == 1 )   CALL dyn_spg_cla( kt )      ! Cross Land Advection (Update (ua,va))
220#endif
221
222      !                                                ! ===============
223      DO jj = 2, jpjm1                                 !  Vertical slab
224         !                                             ! ===============
225
226         ! 2. compute the next vertically averaged velocity
227         ! ------------------------------------------------
228         !     (effect of the additional force not included)
229         ! initialize to zero
230         DO ji = 2, jpim1
231            spgu(ji,jj) = 0.e0
232            spgv(ji,jj) = 0.e0
233         END DO
234
235         ! vertical sum
236         DO jk = 1, jpk
237            DO ji = 2, jpim1
238               spgu(ji,jj) = spgu(ji,jj) + fse3u(ji,jj,jk) * ua(ji,jj,jk)
239               spgv(ji,jj) = spgv(ji,jj) + fse3v(ji,jj,jk) * va(ji,jj,jk)
240            END DO
241         END DO
242
243         ! transport: multiplied by the horizontal scale factor
244         DO ji = 2, jpim1
245            spgu(ji,jj) = spgu(ji,jj) * e2u(ji,jj)
246            spgv(ji,jj) = spgv(ji,jj) * e1v(ji,jj)
247         END DO
248
249         !                                             ! ===============
250      END DO                                           !   End of slab
251      !                                                ! ===============
252
253      !,,,,,,,,,,,,,,,,,,,,,,,,,,,,,synchro,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,
254     
255      ! Boundary conditions on (spgu,spgv)
256      CALL lbc_lnk( spgu, 'U', -1. )
257      CALL lbc_lnk( spgv, 'V', -1. )
258
259      !,,,,,,,,,,,,,,,,,,,,,,,,,,,,,synchro,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,
260
261      ! 3. Right hand side of the elliptic equation and first guess
262      ! -----------------------------------------------------------
263      DO jj = 2, jpjm1
264         DO ji = 2, jpim1
265            ! Divergence of the after vertically averaged velocity
266            zgcb =  spgu(ji,jj) - spgu(ji-1,jj)   &
267                  + spgv(ji,jj) - spgv(ji,jj-1)
268            gcb(ji,jj) = gcdprc(ji,jj) * zgcb
269            ! First guess of the after barotropic transport divergence
270            zbtd = gcx(ji,jj)
271            gcx (ji,jj) = 2. * zbtd   - gcxb(ji,jj)
272            gcxb(ji,jj) =      zbtd
273         END DO
274      END DO
275
276      !,,,,,,,,,,,,,,,,,,,,,,,,,,,,,synchro,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,
277     
278      ! 4. Relative precision (computation on one processor)
279      ! ---------------------
280      rnorme =0.
281      DO jj = 1, jpj
282         DO ji = 1, jpi
283            zgwgt  = gcdmat(ji,jj) * gcb(ji,jj)
284            rnorme = rnorme + gcb(ji,jj) * zgwgt
285         END DO
286      END DO
287      IF( lk_mpp )   CALL mpp_sum( rnorme )   ! sum over the global domain
288
289      epsr = eps * eps * rnorme
290      ncut = 0
291      ! if rnorme is 0, the solution is 0, the solver isn't called
292      IF( rnorme == 0.e0 ) THEN
293         gcx(:,:) = 0.e0
294         res   = 0.e0
295         niter = 0
296         ncut  = 999
297      ENDIF
298     
299      !,,,,,,,,,,,,,,,,,,,,,,,,,,,,,synchro,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,
300     
301      ! 5. Evaluate the next transport divergence
302      ! -----------------------------------------
303      !    Iterarive solver for the elliptic equation (except IF sol.=0)
304      !    (output in gcx with boundary conditions applied)
305      kindic = 0
306      IF( ncut == 0 ) THEN
307         IF( nsolv == 1 ) THEN         ! diagonal preconditioned conjuguate gradient
308            CALL sol_pcg( kindic )
309         ELSEIF( nsolv == 2 ) THEN     ! successive-over-relaxation
310            CALL sol_sor( kindic )
311         ELSEIF( nsolv == 3 ) THEN     ! FETI solver
312            CALL sol_fet( kindic )
313         ELSE                          ! e r r o r in nsolv namelist parameter
314            IF(lwp) WRITE(numout,cform_err)
315            IF(lwp) WRITE(numout,*) ' dyn_spg_fsc_atsk : e r r o r, nsolv = 1, 2 or 3'
316            IF(lwp) WRITE(numout,*) ' ~~~~~~~~~~~~~~~~                not = ', nsolv
317            nstop = nstop + 1
318         ENDIF
319      ENDIF
320
321      !,,,,,,,,,,,,,,,,,,,,,,,,,,,,,synchro,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,
322
323      !                                                ! ===============
324      DO jj = 2, jpjm1                                 !  Vertical slab
325         !                                             ! ===============
326         
327         ! 6. Transport divergence gradient multiplied by z2dt
328         ! -----------------------------------------------====
329         DO ji = 2, jpim1
330            ! trend of Transport divergence gradient
331            ztdgu = znugdt * (gcx(ji+1,jj  ) - gcx(ji,jj) ) / e1u(ji,jj)
332            ztdgv = znugdt * (gcx(ji  ,jj+1) - gcx(ji,jj) ) / e2v(ji,jj)
333            ! multiplied by z2dt
334#if defined key_obc
335            ! caution : grad D = 0 along open boundaries
336            spgu(ji,jj) = z2dt * ztdgu * obcumask(ji,jj)
337            spgv(ji,jj) = z2dt * ztdgv * obcvmask(ji,jj)
338#else
339            spgu(ji,jj) = z2dt * ztdgu
340            spgv(ji,jj) = z2dt * ztdgv
341#endif
342         END DO
343         
344         ! 7.  Add the trends multiplied by z2dt to the after velocity
345         ! -----------------------------------------------------------
346         !     ( c a u t i o n : (ua,va) here are the after velocity not the
347         !                       trend, the leap-frog time stepping will not
348         !                       be done in dynnxt.F routine)
349         DO jk = 1, jpkm1
350            DO ji = 2, jpim1
351               ua(ji,jj,jk) = (ua(ji,jj,jk) + spgu(ji,jj)) * umask(ji,jj,jk)
352               va(ji,jj,jk) = (va(ji,jj,jk) + spgv(ji,jj)) * vmask(ji,jj,jk)
353            END DO
354         END DO
355
356         IF(l_ctl) THEN         ! print sum trends (used for debugging)
357            WRITE(numout,*) ' spg  - Ua: ', SUM( ua(2:nictl,2:njctl,1:jpkm1)*umask(2:nictl,2:njctl,1:jpkm1) ),   &
358               &                   ' Va: ', SUM( va(2:nictl,2:njctl,1:jpkm1)*vmask(2:nictl,2:njctl,1:jpkm1) )
359         ENDIF
360
361         ! 8. Sea surface elevation time stepping
362         ! --------------------------------------
363         ! Euler (forward) time stepping, no time filter
364         IF( neuler == 0 .AND. kt == nit000 ) THEN
365            DO ji = 1, jpi
366               ! after free surface elevation
367               zssha = sshb(ji,jj) + rdt * ( wn(ji,jj,1) - emp(ji,jj) * zraur ) * tmask(ji,jj,1)
368               ! swap of arrays
369               sshb(ji,jj) = sshn(ji,jj)
370               sshn(ji,jj) = zssha
371            END DO
372         ELSE
373            ! Leap-frog time stepping and time filter
374            DO ji = 1, jpi
375               ! after free surface elevation
376               zssha = sshb(ji,jj) + z2dt * ( wn(ji,jj,1) - emp(ji,jj) * zraur ) * tmask(ji,jj,1)
377               ! time filter and array swap
378               sshb(ji,jj) = atfp * ( sshb(ji,jj) + zssha ) + atfp1 * sshn(ji,jj)
379               sshn(ji,jj) = zssha
380            END DO
381         ENDIF
382         !                                             ! ===============
383      END DO                                           !   End of slab
384      !                                                ! ===============
385
386      ! save the surface pressure gradient trends for diagnostic
387      ! momentum trends
388      IF( l_trddyn )   THEN
389         DO jk = 1, jpkm1
390            ztdua(:,:,jk) = ztdua(:,:,jk) + spgu(:,:)/z2dt 
391            ztdva(:,:,jk) = ztdva(:,:,jk) + spgv(:,:)/z2dt 
392         END DO
393
394         CALL trd_mod(ztdua, ztdva, jpdtdspg, 'DYN', kt)
395      ENDIF
396     
397      !Boundary conditions on sshn
398      CALL lbc_lnk( sshn, 'T', 1. )
399
400      IF(l_ctl) THEN         ! print sum trends (used for debugging)
401         WRITE(numout,*) ' spg - ssh:', SUM( sshn(2:nictl,2:njctl)*tmask(2:nictl,2:njctl,1) )
402      ENDIF
403
404
405   END SUBROUTINE dyn_spg_fsc_atsk
406
407#else
408   !!----------------------------------------------------------------------
409   !!   Default case :                                         Empty module
410   !!----------------------------------------------------------------------
411   LOGICAL, PUBLIC, PARAMETER ::   lk_dynspg_fsc_tsk = .FALSE.   !: free surf. cst vol. flag
412CONTAINS
413   SUBROUTINE dyn_spg_fsc_atsk( kt, kindic )      ! Empty module
414      WRITE(*,*) 'dyn_spg_fsc_atsk: You should not have seen this print! error?', kt, kindic
415   END SUBROUTINE dyn_spg_fsc_atsk
416#endif
417   
418   !!======================================================================
419END MODULE dynspg_fsc_atsk
Note: See TracBrowser for help on using the repository browser.