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

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

CT : BUGFIX012 : Running problem for EEL5 configuration is solved

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 17.7 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 trdtra_oce     ! ocean active tracer trend
21   USE trddyn_oce     ! ocean dynamics trend
22   USE in_out_manager  ! I/O manager
23   USE phycst          ! physical constant
24   USE ocesbc          ! Ocen Surface Boundary condition
25   USE flxrnf          ! ???
26   USE sol_oce         ! ocean elliptic solver
27   USE solpcg          ! preconditionned conjugate gradient solver
28   USE solsor          ! Successive Over-relaxation solver
29   USE solfet          ! FETI solver
30   USE obc_oce         ! Lateral open boundary condition
31   USE obcdyn          ! open boudary condition
32   USE obcdyn          !    "              "
33   USE obcvol          !    "              "
34   USE lib_mpp         ! ???
35   USE lbclnk          ! ???
36   USE cla_dynspg      ! ???
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 (utrd,vtrd) ('key_diatrends')
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      !!---------------------------------------------------------------------
109      !! * Arguments
110      INTEGER, INTENT( in ) ::   kt          ! ocean time-step index
111      INTEGER, INTENT( out ) ::   kindic     ! solver convergence flag
112                                             ! if the solver doesn t converge
113                                             ! the flag is < 0
114      !! * Local declarations
115      INTEGER  ::   ji, jj, jk               ! dummy loop indices
116      REAL(wp) ::   &              ! temporary scalars
117         z2dt, z2dtg, zraur, znugdt, znurau,   &
118         zssha, zspgu, zspgv, zegu, zegv,   &
119         zgcb, zbtd, ztdgu, ztdgv, zgwgt
120      !!----------------------------------------------------------------------
121
122      IF( kt == nit000 ) THEN
123         IF(lwp) WRITE(numout,*)
124         IF(lwp) WRITE(numout,*) 'dyn_spg_fsc_atsk : surface pressure gradient trend'
125         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~~~~  (free surface constant volume, autotasking case)'
126
127         ! set to zero free surface specific arrays
128         spgu(:,:) = 0.e0      ! surface pressure gradient (i-direction)
129         spgv(:,:) = 0.e0      ! surface pressure gradient (j-direction)
130      ENDIF
131
132      ! 0. Local constant initialization
133      ! --------------------------------
134      ! time step: leap-frog
135      z2dt = 2. * rdt
136      ! time step: Euler if restart from rest
137      IF( neuler == 0 .AND. kt == nit000 ) z2dt = rdt
138      ! coefficients
139      z2dtg  = grav * z2dt
140      zraur  = 1. / rauw
141      znugdt =  rnu * grav * z2dt
142      znurau =  znugdt * zraur
143      IF( lk_mpp ) THEN
144         ! Mpp : export boundary values of to neighboring processors
145         !!bug ???  why only in mpp?  is it really needed???
146         CALL lbc_lnk( ua, 'U' , -1. )
147         CALL lbc_lnk( va, 'V' , -1. )
148      ENDIF
149
150      !                                                ! ===============
151      DO jj = 2, jpjm1                                 !  Vertical slab
152         !                                             ! ===============
153         ! 1. Surface pressure gradient (now)
154         ! ----------------------------
155         DO ji = 2, jpim1
156            zspgu =   - grav * ( sshn(ji+1,jj) - sshn(ji,jj) ) / e1u(ji,jj)
157            zspgv =   - grav * ( sshn(ji,jj+1) - sshn(ji,jj) ) / e2v(ji,jj)
158            zegu  = + znurau * ( emp (ji+1,jj) - emp (ji,jj) ) / e1u(ji,jj)
159            zegv  = + znurau * ( emp (ji,jj+1) - emp (ji,jj) ) / e2v(ji,jj)
160            spgu(ji,jj) = zspgu + zegu
161            spgv(ji,jj) = zspgv + zegv
162#if defined key_trddyn
163            ! save the surface pressure gradient trends
164            utrd2(ji,jj,1) = zspgu
165            vtrd2(ji,jj,1) = zspgv
166            ! save the first part of the additional force trend
167            utrd2(ji,jj,2) = zegu
168            vtrd2(ji,jj,2) = zegv
169#endif
170         END DO 
171
172         ! 2. Add the surface pressure trend to the general trend
173         ! ------------------------------------------------------
174         DO jk = 1, jpkm1
175            DO ji = 2, jpim1
176               ua(ji,jj,jk) = ua(ji,jj,jk) + spgu(ji,jj)
177               va(ji,jj,jk) = va(ji,jj,jk) + spgv(ji,jj)
178#if defined key_trddyn
179               ! save the surface pressure gradient trend for diagnostics
180               utrd(ji,jj,jk,8) = spgu(ji,jj)
181               vtrd(ji,jj,jk,8) = spgv(ji,jj) 
182#endif
183            END DO
184         END DO
185         !                                             ! ===============
186      END DO                                           !   End of slab
187      !                                                ! ===============
188
189      !,,,,,,,,,,,,,,,,,,,,,,,,,,,,,synchro,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,
190
191      !                                                ! ===============
192      DO jj = 2, jpjm1                                 !  Vertical slab
193         !                                             ! ===============
194
195         ! 1. Evaluate the masked next velocity
196         ! ------------------------------------
197         !     (effect of the additional force not included)
198
199         DO jk = 1, jpkm1
200            DO ji = 2, jpim1
201               ua(ji,jj,jk) = ( ub(ji,jj,jk) + z2dt * ua(ji,jj,jk) ) * umask(ji,jj,jk)
202               va(ji,jj,jk) = ( vb(ji,jj,jk) + z2dt * va(ji,jj,jk) ) * vmask(ji,jj,jk)
203            END DO
204         END DO
205
206         !                                             ! ===============
207      END DO                                           !   End of slab
208      !                                                ! ===============
209
210#if defined key_obc
211         ! Update velocities on each open boundary with the radiation algorithm
212         CALL obc_dyn( kt )
213         ! Correction of the barotropic componant velocity to control the volume of the system
214         CALL obc_vol( kt )
215#endif
216#if defined key_orca_r2
217      IF( n_cla == 1 )   CALL dyn_spg_cla( kt )      ! Cross Land Advection (Update (ua,va))
218#endif
219
220      !                                                ! ===============
221      DO jj = 2, jpjm1                                 !  Vertical slab
222         !                                             ! ===============
223
224         ! 2. compute the next vertically averaged velocity
225         ! ------------------------------------------------
226         !     (effect of the additional force not included)
227         ! initialize to zero
228         DO ji = 2, jpim1
229            spgu(ji,jj) = 0.e0
230            spgv(ji,jj) = 0.e0
231         END DO
232
233         ! vertical sum
234         DO jk = 1, jpk
235            DO ji = 2, jpim1
236               spgu(ji,jj) = spgu(ji,jj) + fse3u(ji,jj,jk) * ua(ji,jj,jk)
237               spgv(ji,jj) = spgv(ji,jj) + fse3v(ji,jj,jk) * va(ji,jj,jk)
238            END DO
239         END DO
240
241         ! transport: multiplied by the horizontal scale factor
242         DO ji = 2, jpim1
243            spgu(ji,jj) = spgu(ji,jj) * e2u(ji,jj)
244            spgv(ji,jj) = spgv(ji,jj) * e1v(ji,jj)
245         END DO
246
247         !                                             ! ===============
248      END DO                                           !   End of slab
249      !                                                ! ===============
250
251      !,,,,,,,,,,,,,,,,,,,,,,,,,,,,,synchro,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,
252     
253      ! Boundary conditions on (spgu,spgv)
254      CALL lbc_lnk( spgu, 'U', -1. )
255      CALL lbc_lnk( spgv, 'V', -1. )
256
257      !,,,,,,,,,,,,,,,,,,,,,,,,,,,,,synchro,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,
258
259      ! 3. Right hand side of the elliptic equation and first guess
260      ! -----------------------------------------------------------
261      DO jj = 2, jpjm1
262         DO ji = 2, jpim1
263            ! Divergence of the after vertically averaged velocity
264            zgcb =  spgu(ji,jj) - spgu(ji-1,jj)   &
265                  + spgv(ji,jj) - spgv(ji,jj-1)
266            gcb(ji,jj) = gcdprc(ji,jj) * zgcb
267            ! First guess of the after barotropic transport divergence
268            zbtd = gcx(ji,jj)
269            gcx (ji,jj) = 2. * zbtd   - gcxb(ji,jj)
270            gcxb(ji,jj) =      zbtd
271         END DO
272      END DO
273
274      !,,,,,,,,,,,,,,,,,,,,,,,,,,,,,synchro,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,
275     
276      ! 4. Relative precision (computation on one processor)
277      ! ---------------------
278      rnorme =0.
279      DO jj = 1, jpj
280         DO ji = 1, jpi
281            zgwgt  = gcdmat(ji,jj) * gcb(ji,jj)
282            rnorme = rnorme + gcb(ji,jj) * zgwgt
283         END DO
284      END DO
285      IF( lk_mpp )   CALL mpp_sum( rnorme )   ! sum over the global domain
286
287      epsr = eps * eps * rnorme
288      ncut = 0
289      ! if rnorme is 0, the solution is 0, the solver isn't called
290      IF( rnorme == 0.e0 ) THEN
291         gcx(:,:) = 0.e0
292         res   = 0.e0
293         niter = 0
294         ncut  = 999
295      ENDIF
296     
297      !,,,,,,,,,,,,,,,,,,,,,,,,,,,,,synchro,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,
298     
299      ! 5. Evaluate the next transport divergence
300      ! -----------------------------------------
301      !    Iterarive solver for the elliptic equation (except IF sol.=0)
302      !    (output in gcx with boundary conditions applied)
303      kindic = 0
304      IF( ncut == 0 ) THEN
305         IF( nsolv == 1 ) THEN         ! diagonal preconditioned conjuguate gradient
306            CALL sol_pcg( kindic )
307         ELSEIF( nsolv == 2 ) THEN     ! successive-over-relaxation
308            CALL sol_sor( kindic )
309         ELSEIF( nsolv == 3 ) THEN     ! FETI solver
310            CALL sol_fet( kindic )
311         ELSE                          ! e r r o r in nsolv namelist parameter
312            IF(lwp) WRITE(numout,cform_err)
313            IF(lwp) WRITE(numout,*) ' dyn_spg_fsc_atsk : e r r o r, nsolv = 1, 2 or 3'
314            IF(lwp) WRITE(numout,*) ' ~~~~~~~~~~~~~~~~                not = ', nsolv
315            nstop = nstop + 1
316         ENDIF
317      ENDIF
318
319      !,,,,,,,,,,,,,,,,,,,,,,,,,,,,,synchro,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,
320
321      !                                                ! ===============
322      DO jj = 2, jpjm1                                 !  Vertical slab
323         !                                             ! ===============
324         
325         ! 6. Transport divergence gradient multiplied by z2dt
326         ! -----------------------------------------------====
327         DO ji = 2, jpim1
328            ! trend of Transport divergence gradient
329            ztdgu = znugdt * (gcx(ji+1,jj  ) - gcx(ji,jj) ) / e1u(ji,jj)
330            ztdgv = znugdt * (gcx(ji  ,jj+1) - gcx(ji,jj) ) / e2v(ji,jj)
331            ! multiplied by z2dt
332#if defined key_obc
333            ! caution : grad D = 0 along open boundaries
334            spgu(ji,jj) = z2dt * ztdgu * obcumask(ji,jj)
335            spgv(ji,jj) = z2dt * ztdgv * obcvmask(ji,jj)
336#else
337            spgu(ji,jj) = z2dt * ztdgu
338            spgv(ji,jj) = z2dt * ztdgv
339#endif
340#if defined key_trddyn
341            ! save the transport divergence gradient trends
342            utrd2(ji,jj,2) = utrd2(ji,jj,2) + ztdgu
343            vtrd2(ji,jj,2) = vtrd2(ji,jj,2) + ztdgv
344#endif
345         END DO
346         
347         ! 7.  Add the trends multiplied by z2dt to the after velocity
348         ! -----------------------------------------------------------
349         !     ( c a u t i o n : (ua,va) here are the after velocity not the
350         !                       trend, the leap-frog time stepping will not
351         !                       be done in dynnxt.F routine)
352         DO jk = 1, jpkm1
353            DO ji = 2, jpim1
354               ua(ji,jj,jk) = (ua(ji,jj,jk) + spgu(ji,jj)) * umask(ji,jj,jk)
355               va(ji,jj,jk) = (va(ji,jj,jk) + spgv(ji,jj)) * vmask(ji,jj,jk)
356#if defined key_trddyn
357               ! save the surface pressure gradient trend for diagnostics
358               utrd(ji,jj,jk,8) = utrd(ji,jj,jk,8) + spgu(ji,jj)/z2dt
359               vtrd(ji,jj,jk,8) = vtrd(ji,jj,jk,8) + spgv(ji,jj)/z2dt
360#endif
361            END DO
362         END DO
363
364         ! 8. Sea surface elevation time stepping
365         ! --------------------------------------
366         ! Euler (forward) time stepping, no time filter
367         IF( neuler == 0 .AND. kt == nit000 ) THEN
368            DO ji = 1, jpi
369               ! after free surface elevation
370               zssha = sshb(ji,jj) + rdt * ( wn(ji,jj,1) - emp(ji,jj) * zraur ) * tmask(ji,jj,1)
371               ! swap of arrays
372               sshb(ji,jj) = sshn(ji,jj)
373               sshn(ji,jj) = zssha
374            END DO
375         ELSE
376            ! Leap-frog time stepping and time filter
377            DO ji = 1, jpi
378               ! after free surface elevation
379               zssha = sshb(ji,jj) + z2dt * ( wn(ji,jj,1) - emp(ji,jj) * zraur ) * tmask(ji,jj,1)
380               ! time filter and array swap
381               sshb(ji,jj) = atfp * ( sshb(ji,jj) + zssha ) + atfp1 * sshn(ji,jj)
382               sshn(ji,jj) = zssha
383            END DO
384         ENDIF
385         !                                             ! ===============
386      END DO                                           !   End of slab
387      !                                                ! ===============
388
389      !Boundary conditions on sshn
390      CALL lbc_lnk( sshn, 'T', 1. )
391
392   END SUBROUTINE dyn_spg_fsc_atsk
393
394#else
395   !!----------------------------------------------------------------------
396   !!   Default case :                                         Empty module
397   !!----------------------------------------------------------------------
398   LOGICAL, PUBLIC, PARAMETER ::   lk_dynspg_fsc_tsk = .FALSE.   !: free surf. cst vol. flag
399CONTAINS
400   SUBROUTINE dyn_spg_fsc_atsk( kt, kindic )      ! Empty module
401      WRITE(*,*) 'dyn_spg_fsc_atsk: You should not have seen this print! error?', kt, kindic
402   END SUBROUTINE dyn_spg_fsc_atsk
403#endif
404   
405   !!======================================================================
406END MODULE dynspg_fsc_atsk
Note: See TracBrowser for help on using the repository browser.