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

source: trunk/NEMO/OPA_SRC/TRA/trazdf_imp_jki.F90 @ 719

Last change on this file since 719 was 719, checked in by ctlod, 17 years ago

get back to the nemo_v2_3 version for trunk

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 10.5 KB
Line 
1MODULE trazdf_imp_jki
2   !!==============================================================================
3   !!                    ***  MODULE  trazdf_imp_jki  ***
4   !! Ocean active tracers:  vertical component of the tracer mixing trend
5   !!==============================================================================
6   !! History :
7   !!   7.0  !  91-11  (G. Madec)  Original code
8   !!        !  92-06  (M. Imbard)  correction on tracer trend loops
9   !!        !  96-01  (G. Madec)  statement function for e3
10   !!        !  97-05  (G. Madec)  vertical component of isopycnal
11   !!        !  97-07  (G. Madec)  geopotential diffusion in s-coord
12   !!        !  00-08  (G. Madec)  double diffusive mixing
13   !!   8.5  !  02-08  (G. Madec)  F90: Free form and module
14   !!   9.0  !  04-08  (C. Talandier) New trends organization
15   !!        !  05-11  (G. Madec)  New step reorganisation
16   !!----------------------------------------------------------------------
17   !!   tra_zdf_imp_jki  : Update the tracer trend with the diagonal vertical
18   !!                 part of the mixing tensor.
19   !!                 use j-k-i loops.
20   !!----------------------------------------------------------------------
21   !! * Modules used
22   USE oce             ! ocean dynamics and tracers variables
23   USE dom_oce         ! ocean space and time domain variables
24   USE ldfslp          ! Make iso-neutral slopes available
25   USE ldftra_oce      ! ocean active tracers: lateral physics
26   USE zdf_oce         ! ocean vertical physics
27   USE zdfddm          ! ocean vertical physics: double diffusion
28   USE trdmod          ! ocean active tracers trends
29   USE trdmod_oce      ! ocean variables trends
30   USE in_out_manager  ! I/O manager
31   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
32   USE prtctl          ! Print control
33
34   IMPLICIT NONE
35   PRIVATE
36
37   !! * Accessibility
38   PUBLIC tra_zdf_imp_jki    ! routine called by step.F90
39
40   !! * Substitutions
41#  include "domzgr_substitute.h90"
42#  include "ldftra_substitute.h90"
43#  include "zdfddm_substitute.h90"
44   !!----------------------------------------------------------------------
45   !!  OPA 9.0 , LOCEAN-IPSL (2005)
46   !!----------------------------------------------------------------------
47
48CONTAINS
49
50   SUBROUTINE tra_zdf_imp_jki( kt, p2dt )
51      !!----------------------------------------------------------------------
52      !!                  ***  ROUTINE tra_zdf_imp_jki  ***
53      !!
54      !! ** Purpose :
55      !!     Compute the trend due to the vertical tracer diffusion inclu-
56      !!     ding the vertical component of lateral mixing (only for second
57      !!     order operator, for fourth order it is already computed and
58      !!     add to the general trend in traldf.F) and add it to the general
59      !!     trend of the tracer equations.
60      !!
61      !! ** Method :
62      !!      save avt coef. resulting from vertical physics alone in zavt:
63      !!         zavt = avt
64      !!      update and save in zavt the vertical eddy viscosity coefficient:
65      !!         avt = avt + wslpi^2+wslj^2
66      !!
67      !!      Second part: vertical trend associated with the vertical physics
68      !!      ===========  (including the vertical flux proportional to dk[t]
69      !!                  associated with the lateral mixing, through the
70      !!                  update of avt)
71      !!      The vertical diffusion of tracers (t & s) is given by:
72      !!             difft = dz( avt dz(t) ) = 1/e3t dk+1( avt/e3w dk(t) )
73      !!      It is computed using a backward time scheme, t=ta.
74      !!      Surface and bottom boundary conditions: no diffusive flux on
75      !!      both tracers (bottom, applied through the masked field avt).
76      !!      Add this trend to the general trend ta,sa :
77      !!         ta = ta + dz( avt dz(t) )
78      !!         (sa = sa + dz( avs dz(t) ) if lk_zdfddm=T )
79      !!
80      !!      Third part: recover avt resulting from the vertical physics
81      !!      ==========  alone, for further diagnostics (for example to
82      !!                  compute the turbocline depth in zdfmxl.F90).
83      !!         avt = zavt
84      !!         (avs = zavs if lk_zdfddm=T )
85      !!
86      !! ** Action :
87      !!         Update (ta,sa) arrays with the before vertical diffusion trend
88      !!
89      !!---------------------------------------------------------------------
90      !! * Modules used
91      USE oce                ,   &
92# if defined key_zdfddm
93         zavs => va,  &  ! use va as workspace
94# endif
95         zavt => ua      ! use ua as workspace
96
97      !! * Arguments
98      INTEGER, INTENT( in ) ::   kt           ! ocean time-step index
99      REAL(wp), DIMENSION(jpk), INTENT( in ) ::   &
100         p2dt                                 ! vertical profile of tracer time-step
101
102      !! * Local declarations
103      INTEGER  ::   ji, jj, jk                ! dummy loop indices
104      INTEGER  ::   ikst, ikenm2, ikstp1      ! temporary integers
105      REAL(wp) ::   zavi, zavip1              ! temporary scalar
106      REAL(wp), DIMENSION(jpi,jpk) ::   &
107         zwd, zws, zwi,         &  ! ???
108         zwx, zwy, zwz, zwt        ! ???
109      !!---------------------------------------------------------------------
110
111      IF( kt == nit000 ) THEN
112         IF(lwp) WRITE(numout,*)
113         IF(lwp) WRITE(numout,*) 'tra_zdf_imp_jki : implicit vertical mixing (j-k-i loops)'
114         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~~'
115      ENDIF
116
117!CDIR PARALLEL DO PRIVATE( zwd, zws, zwi, zwx, zwy, zwz, zwt )
118!$OMP PARALLEL DO PRIVATE( zwd, zws, zwi, zwx, zwy, zwz, zwt )
119      !                                                ! ===============
120      DO jj = 2, jpjm1                                 !  Vertical slab
121         !                                             ! ===============
122
123         ! II. Vertical trend associated with the vertical physics
124         ! =======================================================
125         !     (including the vertical flux proportional to dk[t] associated
126         !      with the lateral mixing, through the avt update)
127         !     dk[ avt dk[ (t,s) ] ] diffusive trends
128
129
130         ! II.0 Matrix construction
131         ! ------------------------
132
133         ! Diagonal, inferior, superior
134         ! (including the bottom boundary condition via avt masked)
135         DO jk = 1, jpkm1
136            DO ji = 2, jpim1
137#if defined key_ldfslp
138               zavi   = avt(ji,jj,jk  ) + fsahtw(ji,jj,jk  )*( wslpi(ji,jj,jk  )*wslpi(ji,jj,jk  )   &
139                  &                                           +wslpj(ji,jj,jk  )*wslpj(ji,jj,jk  ) )
140               zavip1 = avt(ji,jj,jk+1) + fsahtw(ji,jj,jk+1)*( wslpi(ji,jj,jk+1)*wslpi(ji,jj,jk+1)   &
141                  &                                           +wslpj(ji,jj,jk+1)*wslpj(ji,jj,jk+1) )
142#else
143               zavi   = avt(ji,jj,jk  )
144               zavip1 = avt(ji,jj,jk+1)
145#endif
146               zwi(ji,jk) = - p2dt(jk) * zavi   / ( fse3t(ji,jj,jk) * fse3w(ji,jj,jk  ) )
147               zws(ji,jk) = - p2dt(jk) * zavip1 / ( fse3t(ji,jj,jk) * fse3w(ji,jj,jk+1) )
148               zwd(ji,jk) = 1. - zwi(ji,jk) - zws(ji,jk)
149            END DO
150         END DO
151
152         ! Surface boudary conditions
153         DO ji = 2, jpim1
154            zwi(ji,1) = 0.e0
155            zwd(ji,1) = 1. - zws(ji,1)
156         END DO
157
158
159         ! II.1. Vertical diffusion on t
160         ! ---------------------------
161
162         ! Second member construction
163         DO jk = 1, jpkm1
164            DO ji = 2, jpim1
165               zwy(ji,jk) = tb(ji,jj,jk) + p2dt(jk) * ta(ji,jj,jk)
166            END DO
167         END DO
168
169         ! Matrix inversion from the first level
170         ikst = 1
171#   include "zdf.matrixsolver.h90"
172
173         ! Save the masked temperature after in ta
174         ! (c a u t i o n:  temperature not its trend, Leap-frog scheme done
175         !                  it will not be done in tranxt)
176         DO jk = 1, jpkm1
177            DO ji = 2, jpim1
178               ta(ji,jj,jk) = zwx(ji,jk) * tmask(ji,jj,jk)
179            END DO
180         END DO
181
182
183         ! II.2 Vertical diffusion on s
184         ! ---------------------------
185
186#if defined key_zdfddm
187         ! Rebuild the Matrix as avt /= avs
188
189         ! Diagonal, inferior, superior
190         ! (including the bottom boundary condition via avs masked)
191         DO jk = 1, jpkm1
192            DO ji = 2, jpim1
193#if defined key_ldfslp
194               zavi   = fsavs(ji,jj,jk  ) + fsahtw(ji,jj,jk  )*( wslpi(ji,jj,jk  )*wslpi(ji,jj,jk  )   &
195                  &                                             +wslpj(ji,jj,jk  )*wslpj(ji,jj,jk  ) )
196               zavip1 = fsavs(ji,jj,jk+1) + fsahtw(ji,jj,jk+1)*( wslpi(ji,jj,jk+1)*wslpi(ji,jj,jk+1)   &
197                  &                                             +wslpj(ji,jj,jk+1)*wslpj(ji,jj,jk+1) )
198#else
199               zavi   = fsavs(ji,jj,jk  )
200               zavip1 = fsavs(ji,jj,jk+1)
201#endif
202               zwi(ji,jk) = - p2dt(jk) * zavi   / ( fse3t(ji,jj,jk) * fse3w(ji,jj,jk  ) )
203               zws(ji,jk) = - p2dt(jk) * zavip1 / ( fse3t(ji,jj,jk) * fse3w(ji,jj,jk+1) )
204               zwd(ji,jk) = 1. - zwi(ji,jk) - zws(ji,jk)
205            END DO
206         END DO
207
208         ! Surface boudary conditions
209         DO ji = 2, jpim1
210            zwi(ji,1) = 0.e0
211            zwd(ji,1) = 1. - zws(ji,1)
212         END DO
213#endif
214         ! Second member construction
215         DO jk = 1, jpkm1
216            DO ji = 2, jpim1
217               zwy(ji,jk) = sb(ji,jj,jk) + p2dt(jk) * sa(ji,jj,jk)
218            END DO
219         END DO
220
221         ! Matrix inversion from the first level
222         ikst = 1
223#   include "zdf.matrixsolver.h90"
224
225         ! Save the masked salinity after in sa
226         ! (c a u t i o n:  salinity not its trend, Leap-frog scheme done
227         !                  it will not be done in tranxt)
228         DO jk = 1, jpkm1
229            DO ji = 2, jpim1
230               sa(ji,jj,jk) = zwx(ji,jk)  * tmask(ji,jj,jk)
231            END DO
232         END DO
233
234
235#if defined key_ldfslp
236         ! III. recover the avt (avs) resulting from vertical physics only
237         ! ===============================================================
238
239         DO jk = 2, jpkm1
240            DO ji = 2, jpim1
241               avt(ji,jj,jk) = zavt(ji,jj,jk)
242# if defined key_zdfddm
243               fsavs(ji,jj,jk) = zavs(ji,jj,jk)
244# endif
245            END DO
246         END DO
247#endif
248         !                                             ! ===============
249      END DO                                           !   End of slab
250      !                                                ! ===============
251
252   END SUBROUTINE tra_zdf_imp_jki
253
254   !!==============================================================================
255END MODULE trazdf_imp_jki
Note: See TracBrowser for help on using the repository browser.