1 | MODULE trdmod |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE trdmod *** |
---|
4 | !! Ocean diagnostics: ocean tracers and dynamic trends |
---|
5 | !!===================================================================== |
---|
6 | !! History : 1.0 ! 2004-08 (C. Talandier) Original code |
---|
7 | !! - ! 2005-04 (C. Deltel) Add Asselin trend in the ML budget |
---|
8 | !! 3.3 ! 2010-10 (C. Ethe, G. Madec) reorganisation of initialisation phase |
---|
9 | !! 3.5 ! 2012-02 (G. Madec) add 3D trends output for T, S, U, V, PE and KE |
---|
10 | !!---------------------------------------------------------------------- |
---|
11 | #if defined key_trdtra || defined key_trddyn || defined key_trdmld || defined key_trdvor || defined key_esopa |
---|
12 | !!---------------------------------------------------------------------- |
---|
13 | !! trd_mod : manage the type of trend diagnostics |
---|
14 | !! trd_3Diom : output 3D momentum and/or tracer trends using IOM |
---|
15 | !! trd_budget : domain averaged budget of trends (including kinetic energy and tracer variance trends) |
---|
16 | !! trd_mod_init : Initialization step |
---|
17 | !!---------------------------------------------------------------------- |
---|
18 | USE oce ! ocean dynamics and tracers variables |
---|
19 | USE dom_oce ! ocean space and time domain variables |
---|
20 | USE zdf_oce ! ocean vertical physics variables |
---|
21 | USE trdmod_oce ! ocean variables trends |
---|
22 | USE ldftra_oce ! ocean active tracers lateral physics |
---|
23 | USE sbc_oce ! surface boundary condition: ocean |
---|
24 | USE phycst ! physical constants |
---|
25 | USE trdvor ! ocean vorticity trends |
---|
26 | USE trdicp ! ocean bassin integral constraints properties |
---|
27 | USE trdmld ! ocean active mixed layer tracers trends |
---|
28 | USE in_out_manager ! I/O manager |
---|
29 | USE iom ! I/O manager library |
---|
30 | USE lib_mpp ! MPP library |
---|
31 | USE wrk_nemo ! Memory allocation |
---|
32 | |
---|
33 | IMPLICIT NONE |
---|
34 | PRIVATE |
---|
35 | |
---|
36 | REAL(wp) :: r2dt ! time-step, = 2 rdttra except at nit000 (=rdttra) if neuler=0 |
---|
37 | |
---|
38 | PUBLIC trd_mod ! called by all dynXX or traXX modules |
---|
39 | PUBLIC trd_mod_init ! called by opa.F90 module |
---|
40 | |
---|
41 | !! * Substitutions |
---|
42 | # include "domzgr_substitute.h90" |
---|
43 | # include "vectopt_loop_substitute.h90" |
---|
44 | !!---------------------------------------------------------------------- |
---|
45 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
46 | !! $Id$ |
---|
47 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
48 | !!---------------------------------------------------------------------- |
---|
49 | CONTAINS |
---|
50 | |
---|
51 | SUBROUTINE trd_mod( ptrdx, ptrdy, ktrd, ctype, kt ) |
---|
52 | !!--------------------------------------------------------------------- |
---|
53 | !! *** ROUTINE trd_mod *** |
---|
54 | !! |
---|
55 | !! ** Purpose : Dispatch all trends computation, e.g. 3D output, integral |
---|
56 | !! constraints, barotropic vorticity, kinetic enrgy, |
---|
57 | !! potential energy, and/or mixed layer budget. |
---|
58 | !!---------------------------------------------------------------------- |
---|
59 | REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: ptrdx ! Temperature or U trend |
---|
60 | REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: ptrdy ! Salinity or V trend |
---|
61 | INTEGER , INTENT(in ) :: ktrd ! tracer trend index |
---|
62 | CHARACTER(len=3) , INTENT(in ) :: ctype ! momentum or tracers trends type 'DYN'/'TRA' |
---|
63 | INTEGER , INTENT(in ) :: kt ! time step |
---|
64 | !! |
---|
65 | INTEGER :: ji, jj ! dummy loop indices |
---|
66 | REAL(wp), POINTER, DIMENSION(:,:) :: ztswu, ztswv ! 2D workspace |
---|
67 | !!---------------------------------------------------------------------- |
---|
68 | |
---|
69 | CALL wrk_alloc( jpi, jpj, ztswu, ztswv ) |
---|
70 | |
---|
71 | IF( neuler == 0 .AND. kt == nit000 ) THEN ; r2dt = rdt ! = rdtra (restart with Euler time stepping) |
---|
72 | ELSEIF( kt <= nit000 + 1) THEN ; r2dt = 2. * rdt ! = 2 rdttra (leapfrog) |
---|
73 | ENDIF |
---|
74 | |
---|
75 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
76 | IF( ln_3D_trd_d .OR. ln_3D_trd_t ) THEN ! 3D output of momentum and/or tracers trends using IOM interface |
---|
77 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
78 | CALL trd_3Diom ( ptrdx, ptrdy, ktrd, ctype, kt ) |
---|
79 | ! |
---|
80 | ENDIF |
---|
81 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
82 | IF( ln_glo_trd ) THEN ! I. Integral Constraints Properties for momentum and/or tracers trends |
---|
83 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
84 | CALL trd_budget( ptrdx, ptrdy, ktrd, ctype, kt ) |
---|
85 | ! |
---|
86 | ENDIF |
---|
87 | |
---|
88 | ! !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
89 | IF( lk_trdvor .AND. ctype == 'DYN' ) THEN ! II. Vorticity trends |
---|
90 | ! !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
91 | SELECT CASE ( ktrd ) |
---|
92 | CASE ( jpdyn_trd_hpg ) ; CALL trd_vor_zint( ptrdx, ptrdy, jpvor_prg ) ! Hydrostatique Pressure Gradient |
---|
93 | CASE ( jpdyn_trd_keg ) ; CALL trd_vor_zint( ptrdx, ptrdy, jpvor_keg ) ! KE Gradient |
---|
94 | CASE ( jpdyn_trd_rvo ) ; CALL trd_vor_zint( ptrdx, ptrdy, jpvor_rvo ) ! Relative Vorticity |
---|
95 | CASE ( jpdyn_trd_pvo ) ; CALL trd_vor_zint( ptrdx, ptrdy, jpvor_pvo ) ! Planetary Vorticity Term |
---|
96 | CASE ( jpdyn_trd_ldf ) ; CALL trd_vor_zint( ptrdx, ptrdy, jpvor_ldf ) ! Horizontal Diffusion |
---|
97 | CASE ( jpdyn_trd_zad ) ; CALL trd_vor_zint( ptrdx, ptrdy, jpvor_zad ) ! Vertical Advection |
---|
98 | CASE ( jpdyn_trd_spg ) ; CALL trd_vor_zint( ptrdx, ptrdy, jpvor_spg ) ! Surface Pressure Grad. |
---|
99 | CASE ( jpdyn_trd_zdf ) ! Vertical Diffusion |
---|
100 | ztswu(:,:) = 0.e0 ; ztswv(:,:) = 0.e0 |
---|
101 | DO jj = 2, jpjm1 ! wind stress trends |
---|
102 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
103 | ztswu(ji,jj) = ( utau_b(ji,jj) + utau(ji,jj) ) / ( fse3u(ji,jj,1) * rau0 ) |
---|
104 | ztswv(ji,jj) = ( vtau_b(ji,jj) + vtau(ji,jj) ) / ( fse3v(ji,jj,1) * rau0 ) |
---|
105 | END DO |
---|
106 | END DO |
---|
107 | ! |
---|
108 | CALL trd_vor_zint( ptrdx, ptrdy, jpvor_zdf ) ! zdf trend including surf./bot. stresses |
---|
109 | CALL trd_vor_zint( ztswu, ztswv, jpvor_swf ) ! surface wind stress |
---|
110 | CASE ( jpdyn_trd_bfr ) |
---|
111 | CALL trd_vor_zint( ptrdx, ptrdy, jpvor_bfr ) ! Bottom stress |
---|
112 | END SELECT |
---|
113 | ! |
---|
114 | ENDIF |
---|
115 | |
---|
116 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
117 | ! III. Mixed layer trends for active tracers |
---|
118 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
119 | |
---|
120 | IF( lk_trdmld .AND. ctype == 'TRA' ) THEN |
---|
121 | !----------------------------------------------------------------------------------------------- |
---|
122 | ! W.A.R.N.I.N.G : |
---|
123 | ! jptra_trd_ldf : called by traldf.F90 |
---|
124 | ! at this stage we store: |
---|
125 | ! - the lateral geopotential diffusion (here, lateral = horizontal) |
---|
126 | ! - and the iso-neutral diffusion if activated |
---|
127 | ! jptra_trd_zdf : called by trazdf.F90 |
---|
128 | ! * in case of iso-neutral diffusion we store the vertical diffusion component in the |
---|
129 | ! lateral trend including the K_z contrib, which will be removed later (see trd_mld) |
---|
130 | !----------------------------------------------------------------------------------------------- |
---|
131 | |
---|
132 | SELECT CASE ( ktrd ) |
---|
133 | CASE ( jptra_trd_xad ) ; CALL trd_mld_zint( ptrdx, ptrdy, jpmld_xad, '3D' ) ! zonal advection |
---|
134 | CASE ( jptra_trd_yad ) ; CALL trd_mld_zint( ptrdx, ptrdy, jpmld_yad, '3D' ) ! merid. advection |
---|
135 | CASE ( jptra_trd_zad ) ; CALL trd_mld_zint( ptrdx, ptrdy, jpmld_zad, '3D' ) ! vertical advection |
---|
136 | CASE ( jptra_trd_ldf ) ; CALL trd_mld_zint( ptrdx, ptrdy, jpmld_ldf, '3D' ) ! lateral diffusion |
---|
137 | CASE ( jptra_trd_bbl ) ; CALL trd_mld_zint( ptrdx, ptrdy, jpmld_bbl, '3D' ) ! bottom boundary layer |
---|
138 | CASE ( jptra_trd_zdf ) |
---|
139 | IF( ln_traldf_iso ) THEN ; CALL trd_mld_zint( ptrdx, ptrdy, jpmld_ldf, '3D' ) ! lateral diffusion (K_z) |
---|
140 | ELSE ; CALL trd_mld_zint( ptrdx, ptrdy, jpmld_zdf, '3D' ) ! vertical diffusion (K_z) |
---|
141 | ENDIF |
---|
142 | CASE ( jptra_trd_dmp ) ; CALL trd_mld_zint( ptrdx, ptrdy, jpmld_dmp, '3D' ) ! internal 3D restoring (tradmp) |
---|
143 | CASE ( jptra_trd_qsr ) ; CALL trd_mld_zint( ptrdx, ptrdy, jpmld_for, '3D' ) ! air-sea : penetrative sol radiat |
---|
144 | CASE ( jptra_trd_nsr ) |
---|
145 | ptrdx(:,:,2:jpk) = 0._wp ; ptrdy(:,:,2:jpk) = 0._wp |
---|
146 | CALL trd_mld_zint( ptrdx, ptrdy, jpmld_for, '2D' ) ! air-sea : non penetr sol radiat |
---|
147 | CASE ( jptra_trd_bbc ) ; CALL trd_mld_zint( ptrdx, ptrdy, jpmld_bbc, '3D' ) ! bottom bound cond (geoth flux) |
---|
148 | CASE ( jptra_trd_atf ) ; CALL trd_mld_zint( ptrdx, ptrdy, jpmld_atf, '3D' ) ! asselin numerical |
---|
149 | CASE ( jptra_trd_npc ) ; CALL trd_mld_zint( ptrdx, ptrdy, jpmld_npc, '3D' ) ! non penetr convect adjustment |
---|
150 | END SELECT |
---|
151 | ! |
---|
152 | ENDIF |
---|
153 | ! |
---|
154 | CALL wrk_dealloc( jpi, jpj, ztswu, ztswv ) |
---|
155 | ! |
---|
156 | END SUBROUTINE trd_mod |
---|
157 | |
---|
158 | |
---|
159 | SUBROUTINE trd_budget( ptrdx, ptrdy, ktrd, ctype, kt ) |
---|
160 | !!--------------------------------------------------------------------- |
---|
161 | !! *** ROUTINE trd_budget *** |
---|
162 | !! |
---|
163 | !! ** Purpose : integral constraint diagnostics for momentum and/or tracer trends |
---|
164 | !! |
---|
165 | !!---------------------------------------------------------------------- |
---|
166 | REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: ptrdx ! Temperature or U trend |
---|
167 | REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: ptrdy ! Salinity or V trend |
---|
168 | INTEGER , INTENT(in ) :: ktrd ! tracer trend index |
---|
169 | CHARACTER(len=3) , INTENT(in ) :: ctype ! momentum or tracers trends type 'DYN'/'TRA' |
---|
170 | INTEGER , INTENT(in ) :: kt ! time step |
---|
171 | !! |
---|
172 | INTEGER :: ji, jj ! dummy loop indices |
---|
173 | REAL(wp), POINTER, DIMENSION(:,:) :: ztswu, ztswv, z2dx, z2dy ! 2D workspace |
---|
174 | !!---------------------------------------------------------------------- |
---|
175 | |
---|
176 | CALL wrk_alloc( jpi, jpj, ztswu, ztswv, z2dx, z2dy ) |
---|
177 | |
---|
178 | IF( MOD(kt,nn_trd) == 0 .OR. kt == nit000 .OR. kt == nitend ) THEN |
---|
179 | ! |
---|
180 | IF( lk_trdtra .AND. ctype == 'TRA' ) THEN ! active tracer trends |
---|
181 | SELECT CASE ( ktrd ) |
---|
182 | CASE ( jptra_trd_ldf ) ; CALL trd_icp( ptrdx, ptrdy, jpicpt_ldf, ctype ) ! lateral diff |
---|
183 | CASE ( jptra_trd_zdf ) ; CALL trd_icp( ptrdx, ptrdy, jpicpt_zdf, ctype ) ! vertical diff (Kz) |
---|
184 | CASE ( jptra_trd_bbc ) ; CALL trd_icp( ptrdx, ptrdy, jpicpt_bbc, ctype ) ! bottom boundary cond |
---|
185 | CASE ( jptra_trd_bbl ) ; CALL trd_icp( ptrdx, ptrdy, jpicpt_bbl, ctype ) ! bottom boundary layer |
---|
186 | CASE ( jptra_trd_npc ) ; CALL trd_icp( ptrdx, ptrdy, jpicpt_npc, ctype ) ! static instability mixing |
---|
187 | CASE ( jptra_trd_dmp ) ; CALL trd_icp( ptrdx, ptrdy, jpicpt_dmp, ctype ) ! damping |
---|
188 | CASE ( jptra_trd_qsr ) ; CALL trd_icp( ptrdx, ptrdy, jpicpt_qsr, ctype ) ! penetrative solar radiat. |
---|
189 | CASE ( jptra_trd_nsr ) ; z2dx(:,:) = ptrdx(:,:,1) ! non solar radiation |
---|
190 | z2dy(:,:) = ptrdy(:,:,1) |
---|
191 | CALL trd_icp( z2dx , z2dy , jpicpt_nsr, ctype ) |
---|
192 | CASE ( jptra_trd_xad ) ; CALL trd_icp( ptrdx, ptrdy, jpicpt_xad, ctype ) ! x- horiz adv |
---|
193 | CASE ( jptra_trd_yad ) ; CALL trd_icp( ptrdx, ptrdy, jpicpt_yad, ctype ) ! y- horiz adv |
---|
194 | CASE ( jptra_trd_zad ) ; CALL trd_icp( ptrdx, ptrdy, jpicpt_zad, ctype ) ! z- vertical adv |
---|
195 | ! ! surface flux |
---|
196 | IF( lk_vvl ) THEN ! variable volume = zero |
---|
197 | z2dx(:,:) = 0._wp |
---|
198 | z2dy(:,:) = 0._wp |
---|
199 | ELSE ! constant volume = wn*tsn/e3t |
---|
200 | z2dx(:,:) = wn(:,:,1) * tsn(:,:,1,jp_tem) / fse3t(:,:,1) |
---|
201 | z2dy(:,:) = wn(:,:,1) * tsn(:,:,1,jp_sal) / fse3t(:,:,1) |
---|
202 | ENDIF |
---|
203 | CALL trd_icp( z2dx , z2dy , jpicpt_zl1, ctype ) |
---|
204 | END SELECT |
---|
205 | ENDIF |
---|
206 | |
---|
207 | IF( lk_trddyn .AND. ctype == 'DYN' ) THEN ! momentum trends |
---|
208 | ! |
---|
209 | SELECT CASE ( ktrd ) |
---|
210 | CASE( jpdyn_trd_hpg ) ; CALL trd_icp( ptrdx, ptrdy, jpicpd_hpg, ctype ) ! hydrost. pressure gradient |
---|
211 | CASE( jpdyn_trd_spg ) ; CALL trd_icp( ptrdx, ptrdy, jpicpd_spg, ctype ) ! surface pressure grad. |
---|
212 | CASE( jpdyn_trd_pvo ) ; CALL trd_icp( ptrdx, ptrdy, jpicpd_pvo, ctype ) ! planetary vorticity |
---|
213 | CASE( jpdyn_trd_rvo ) ; CALL trd_icp( ptrdx, ptrdy, jpicpd_rvo, ctype ) ! relative vorticity or metric term |
---|
214 | CASE( jpdyn_trd_keg ) ; CALL trd_icp( ptrdx, ptrdy, jpicpd_keg, ctype ) ! KE gradient or hor. advection |
---|
215 | CASE( jpdyn_trd_zad ) ; CALL trd_icp( ptrdx, ptrdy, jpicpd_zad, ctype ) ! vertical advection |
---|
216 | CASE( jpdyn_trd_ldf ) ; CALL trd_icp( ptrdx, ptrdy, jpicpd_ldf, ctype ) ! lateral diffusion |
---|
217 | CASE( jpdyn_trd_zdf ) ; CALL trd_icp( ptrdx, ptrdy, jpicpd_zdf, ctype ) ! vertical diffusion (icluding bfr & tau) |
---|
218 | ztswu(:,:) = ( utau_b(ji,jj) + utau(ji,jj) ) / ( fse3u(:,:,1) * rau0 ) |
---|
219 | ztswv(:,:) = ( vtau_b(ji,jj) + vtau(ji,jj) ) / ( fse3v(:,:,1) * rau0 ) |
---|
220 | CALL trd_icp( ztswu, ztswv, jpicpd_swf, ctype ) ! wind stress trends |
---|
221 | CASE( jpdyn_trd_bfr ) ; CALL trd_icp( ptrdx, ptrdy, jpicpd_bfr, ctype ) ! bottom friction trends |
---|
222 | END SELECT |
---|
223 | ! |
---|
224 | ENDIF |
---|
225 | ! |
---|
226 | ENDIF |
---|
227 | ! |
---|
228 | CALL wrk_dealloc( jpi, jpj, ztswu, ztswv, z2dx, z2dy ) |
---|
229 | ! |
---|
230 | END SUBROUTINE trd_budget |
---|
231 | |
---|
232 | |
---|
233 | SUBROUTINE trd_3Diom( ptrdx, ptrdy, ktrd, ctype, kt ) |
---|
234 | !!--------------------------------------------------------------------- |
---|
235 | !! *** ROUTINE trd_3Diom *** |
---|
236 | !! |
---|
237 | !! ** Purpose : output 3D trends using IOM |
---|
238 | !!---------------------------------------------------------------------- |
---|
239 | REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: ptrdx ! Temperature or U trend |
---|
240 | REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: ptrdy ! Salinity or V trend |
---|
241 | INTEGER , INTENT(in ) :: ktrd ! tracer trend index |
---|
242 | CHARACTER(len=3) , INTENT(in ) :: ctype ! momentum or tracers trends type 'DYN'/'TRA' |
---|
243 | INTEGER , INTENT(in ) :: kt ! time step |
---|
244 | !! |
---|
245 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
246 | REAL(wp), POINTER, DIMENSION(:,:) :: z2dx, z2dy, ztswu, ztswv ! 2D workspace |
---|
247 | REAL(wp), POINTER, DIMENSION(:,:,:) :: z3dx, z3dy ! 3D workspace |
---|
248 | !!---------------------------------------------------------------------- |
---|
249 | |
---|
250 | IF( lk_trdtra .AND. ctype == 'TRA' ) THEN ! active tracer trends |
---|
251 | ! |
---|
252 | !!gm Rq: mask the trends already masked in trd_tra, but lbc_lnk should probably be added |
---|
253 | ! |
---|
254 | SELECT CASE( ktrd ) |
---|
255 | CASE( jptra_trd_xad ) ; CALL iom_put( "ttrd_xad", ptrdx ) ! x- horizontal advection |
---|
256 | CALL iom_put( "strd_xad", ptrdy ) |
---|
257 | CASE( jptra_trd_yad ) ; CALL iom_put( "ttrd_yad", ptrdx ) ! y- horizontal advection |
---|
258 | CALL iom_put( "strd_yad", ptrdy ) |
---|
259 | CASE( jptra_trd_zad ) ; CALL iom_put( "ttrd_zad", ptrdx ) ! z- vertical advection |
---|
260 | CALL iom_put( "strd_zad", ptrdy ) |
---|
261 | IF( .NOT.lk_vvl ) THEN ! cst volume : adv flux through z=0 surface |
---|
262 | z2dx(:,:) = wn(:,:,1) * tsn(:,:,1,jp_tem) / fse3t(:,:,1) |
---|
263 | z2dy(:,:) = wn(:,:,1) * tsn(:,:,1,jp_sal) / fse3t(:,:,1) |
---|
264 | CALL iom_put( "ttrd_sad", z2dx ) |
---|
265 | CALL iom_put( "strd_sad", z2dy ) |
---|
266 | ENDIF |
---|
267 | CASE( jptra_trd_ldf ) ; CALL iom_put( "ttrd_ldf", ptrdx ) ! lateral diffusion |
---|
268 | CALL iom_put( "strd_ldf", ptrdy ) |
---|
269 | CASE( jptra_trd_zdf ) ; CALL iom_put( "ttrd_zdf", ptrdx ) ! vertical diffusion (including Kz contribution) |
---|
270 | CALL iom_put( "strd_zdf", ptrdy ) |
---|
271 | CASE( jptra_trd_dmp ) ; CALL iom_put( "ttrd_dmp", ptrdx ) ! internal restoring (damping) |
---|
272 | CALL iom_put( "strd_dmp", ptrdy ) |
---|
273 | CASE( jptra_trd_bbl ) ; CALL iom_put( "ttrd_bbl", ptrdx ) ! bottom boundary layer |
---|
274 | CALL iom_put( "strd_bbl", ptrdy ) |
---|
275 | CASE( jptra_trd_npc ) ; CALL iom_put( "ttrd_npc", ptrdx ) ! static instability mixing |
---|
276 | CALL iom_put( "strd_npc", ptrdy ) |
---|
277 | CASE( jptra_trd_qsr ) ; CALL iom_put( "ttrd_qsr", ptrdx ) ! penetrative solar radiat. (only on temperature) |
---|
278 | CASE( jptra_trd_nsr ) ; CALL iom_put( "ttrd_qns", ptrdx(:,:,1) ) ! non-solar radiation (only on temperature) |
---|
279 | CASE( jptra_trd_bbc ) ; CALL iom_put( "ttrd_bbc", ptrdx ) ! geothermal heating (only on temperature) |
---|
280 | END SELECT |
---|
281 | ENDIF |
---|
282 | |
---|
283 | IF( lk_trddyn .AND. ctype == 'DYN' ) THEN ! momentum trends |
---|
284 | ! |
---|
285 | ptrdx(:,:,:) = ptrdx(:,:,:) * umask(:,:,:) ! mask the trends |
---|
286 | ptrdy(:,:,:) = ptrdy(:,:,:) * vmask(:,:,:) |
---|
287 | !!gm NB : here a lbc_lnk should probably be added |
---|
288 | ! |
---|
289 | SELECT CASE( ktrd ) |
---|
290 | CASE( jpdyn_trd_hpg ) ; CALL iom_put( "utrd_hpg", ptrdx ) ! hydrostatic pressure gradient |
---|
291 | CALL iom_put( "vtrd_hpg", ptrdy ) |
---|
292 | CASE( jpdyn_trd_spg ) ; CALL iom_put( "utrd_spg", ptrdx ) ! surface pressure gradient |
---|
293 | CALL iom_put( "vtrd_spg", ptrdy ) |
---|
294 | CASE( jpdyn_trd_pvo ) ; CALL iom_put( "utrd_pvo", ptrdx ) ! planetary vorticity |
---|
295 | CALL iom_put( "vtrd_pvo", ptrdy ) |
---|
296 | CASE( jpdyn_trd_rvo ) ; CALL iom_put( "utrd_rvo", ptrdx ) ! relative vorticity (or metric term) |
---|
297 | CALL iom_put( "vtrd_rvo", ptrdy ) |
---|
298 | CASE( jpdyn_trd_keg ) ; CALL iom_put( "utrd_keg", ptrdx ) ! Kinetic Energy gradient (or had) |
---|
299 | CALL iom_put( "vtrd_keg", ptrdy ) |
---|
300 | z3dx(:,:,:) = 0._wp ! U.dxU & V.dyV (approximation) |
---|
301 | z3dy(:,:,:) = 0._wp |
---|
302 | DO jk = 1, jpkm1 ! no mask as un,vn are masked |
---|
303 | DO jj = 2, jpjm1 |
---|
304 | DO ji = 2, jpim1 |
---|
305 | z3dx(ji,jj,jk) = un(ji,jj,jk) * ( un(ji+1,jj,jk) - un(ji-1,jj,jk) ) / ( 2._wp * e1u(ji,jj) ) |
---|
306 | z3dy(ji,jj,jk) = vn(ji,jj,jk) * ( vn(ji,jj+1,jk) - vn(ji,jj-1,jk) ) / ( 2._wp * e2v(ji,jj) ) |
---|
307 | END DO |
---|
308 | END DO |
---|
309 | END DO |
---|
310 | CALL lbc_lnk( z3dx, 'U', -1. ) ; CALL lbc_lnk( z3dy, 'V', -1. ) |
---|
311 | CALL iom_put( "utrd_udx", z3dx ) |
---|
312 | CALL iom_put( "vtrd_vdy", z3dy ) |
---|
313 | CASE( jpdyn_trd_zad ) ; CALL iom_put( "utrd_zad", ptrdx ) ! vertical advection |
---|
314 | CALL iom_put( "vtrd_zad", ptrdy ) |
---|
315 | CASE( jpdyn_trd_ldf ) ; CALL iom_put( "utrd_ldf", ptrdx ) ! lateral diffusion |
---|
316 | CALL iom_put( "vtrd_ldf", ptrdy ) |
---|
317 | CASE( jpdyn_trd_zdf ) ; CALL iom_put( "utrd_zdf", ptrdx ) ! vertical diffusion |
---|
318 | CALL iom_put( "vtrd_zdf", ptrdy ) |
---|
319 | ! ! wind stress trends |
---|
320 | z2dx(:,:) = ( utau_b(ji,jj) + utau(ji,jj) ) / ( fse3u(:,:,1) * rau0 ) |
---|
321 | z2dy(:,:) = ( vtau_b(ji,jj) + vtau(ji,jj) ) / ( fse3v(:,:,1) * rau0 ) |
---|
322 | CALL iom_put( "utrd_tau", z2dx ) |
---|
323 | CALL iom_put( "vtrd_tau", z2dy ) |
---|
324 | CASE( jpdyn_trd_bfr ) ; CALL iom_put( "utrd_bfr", ptrdx ) ! bottom friction term |
---|
325 | CALL iom_put( "vtrd_bfr", ptrdy ) |
---|
326 | END SELECT |
---|
327 | ! |
---|
328 | ENDIF |
---|
329 | ! |
---|
330 | CALL wrk_dealloc( jpi, jpj , z2dx, z2dy, ztswu, ztswv ) |
---|
331 | CALL wrk_dealloc( jpi, jpj, jpk, z3dx, z3dy ) |
---|
332 | ! |
---|
333 | END SUBROUTINE trd_3Diom |
---|
334 | |
---|
335 | #else |
---|
336 | !!---------------------------------------------------------------------- |
---|
337 | !! Default case : Empty module No trend diagnostics |
---|
338 | !!---------------------------------------------------------------------- |
---|
339 | CONTAINS |
---|
340 | SUBROUTINE trd_mod( ptrdx, ptrdy, ktrd, ctype, kt ) ! Empty routine |
---|
341 | REAL :: ptrdx(:,:,:), ptrdy(:,:,:) |
---|
342 | INTEGER :: ktrd, kt |
---|
343 | CHARACTER(len=3) :: ctype |
---|
344 | WRITE(*,*) 'trd_mod: You should not have seen this print! error ?', & |
---|
345 | & ptrdx(1,1,1), ptrdy(1,1,1), ktrd, ctype, kt |
---|
346 | END SUBROUTINE trd_mod |
---|
347 | #endif |
---|
348 | |
---|
349 | SUBROUTINE trd_mod_init |
---|
350 | !!---------------------------------------------------------------------- |
---|
351 | !! *** ROUTINE trd_mod_init *** |
---|
352 | !! |
---|
353 | !! ** Purpose : Initialization of activated trends |
---|
354 | !!---------------------------------------------------------------------- |
---|
355 | USE in_out_manager ! I/O manager |
---|
356 | |
---|
357 | NAMELIST/namtrd/ ln_3D_trd_d, ln_KE_trd, ln_vor_trd, ln_ML_trd_d, & |
---|
358 | & ln_3D_trd_t, ln_PE_trd, ln_glo_trd, ln_ML_trd_t, & |
---|
359 | & nn_trd , cn_trdrst_in , ln_trdmld_restart, & |
---|
360 | & nn_ctls, cn_trdrst_out, ln_trdmld_instant, rn_ucf |
---|
361 | !!---------------------------------------------------------------------- |
---|
362 | |
---|
363 | IF( l_trdtra .OR. l_trddyn ) THEN |
---|
364 | REWIND( numnam ) |
---|
365 | READ ( numnam, namtrd ) ! namelist namtrd : trends diagnostic |
---|
366 | |
---|
367 | IF(lwp) THEN |
---|
368 | WRITE(numout,*) |
---|
369 | WRITE(numout,*) ' trd_mod_init : Momentum/Tracers trends' |
---|
370 | WRITE(numout,*) ' ~~~~~~~~~~~~~' |
---|
371 | WRITE(numout,*) ' Namelist namtrd : set trends parameters' |
---|
372 | WRITE(numout,*) ' U & V trends: 3D output ln_3D_trd_d = ', ln_3D_trd_d |
---|
373 | WRITE(numout,*) ' T & S trends: 3D output ln_3D_trd_t = ', ln_3D_trd_t |
---|
374 | WRITE(numout,*) ' Kinetic Energy trends ln_KE_trd = ', ln_KE_trd |
---|
375 | WRITE(numout,*) ' Potential Energy trends ln_PE_trd = ', ln_PE_trd |
---|
376 | WRITE(numout,*) ' Barotropic vorticity trends ln_vor_trd = ', ln_vor_trd |
---|
377 | WRITE(numout,*) ' check domain averaged dyn & tra trends ln_glo_trd = ', ln_glo_trd |
---|
378 | WRITE(numout,*) ' U & V trends: Mixed Layer averaged ln_ML_trd_d = ', ln_3D_trd_d |
---|
379 | WRITE(numout,*) ' T & S trends: Mixed Layer averaged ln_ML_trd_t = ', ln_3D_trd_t |
---|
380 | ! |
---|
381 | WRITE(numout,*) ' frequency of trends diagnostics (glo) nn_trd = ', nn_trd |
---|
382 | WRITE(numout,*) ' control surface type (mld) nn_ctls = ', nn_ctls |
---|
383 | WRITE(numout,*) ' restart for ML diagnostics ln_trdmld_restart = ', ln_trdmld_restart |
---|
384 | WRITE(numout,*) ' instantaneous or mean ML T/S ln_trdmld_instant = ', ln_trdmld_instant |
---|
385 | WRITE(numout,*) ' unit conversion factor rn_ucf = ', rn_ucf |
---|
386 | ENDIF |
---|
387 | ENDIF |
---|
388 | ! |
---|
389 | IF( ln_KE_trd .OR. ln_PE_trd .OR. ln_ML_trd_d ) & |
---|
390 | CALL ctl_stop( 'KE, PE, aur ML on momentum are not yet coded we stop' ) |
---|
391 | !!gm : Potential BUG : 3D output only for vector invariant form! add a ctl_stop or code the flux form case |
---|
392 | !!gm : bug/pb for vertical advection of tracer in vvl case: add T.dt[eta] in the output... |
---|
393 | ! |
---|
394 | IF( lk_trddyn .OR. lk_trdtra ) CALL trd_icp_init ! integral constraints trends |
---|
395 | IF( lk_trdmld ) CALL trd_mld_init ! mixed-layer trends (active tracers) |
---|
396 | IF( lk_trdvor ) CALL trd_vor_init ! vorticity trends |
---|
397 | ! |
---|
398 | END SUBROUTINE trd_mod_init |
---|
399 | |
---|
400 | !!====================================================================== |
---|
401 | END MODULE trdmod |
---|