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

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

CT : UPDATE151 : New trends organization

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