1 | MODULE traadv |
---|
2 | !!============================================================================== |
---|
3 | !! *** MODULE traadv *** |
---|
4 | !! Ocean active tracers: advection trend |
---|
5 | !!============================================================================== |
---|
6 | !! History : 2.0 ! 2005-11 (G. Madec) Original code |
---|
7 | !! 3.3 ! 2010-09 (C. Ethe, G. Madec) merge TRC-TRA + switch from velocity to transport |
---|
8 | !!---------------------------------------------------------------------- |
---|
9 | |
---|
10 | !!---------------------------------------------------------------------- |
---|
11 | !! tra_adv : compute ocean tracer advection trend |
---|
12 | !! tra_adv_ctl : control the different options of advection scheme |
---|
13 | !!---------------------------------------------------------------------- |
---|
14 | USE oce ! ocean dynamics and active tracers |
---|
15 | USE dom_oce ! ocean space and time domain |
---|
16 | USE traadv_cen2 ! 2nd order centered scheme (tra_adv_cen2 routine) |
---|
17 | USE traadv_tvd ! TVD scheme (tra_adv_tvd routine) |
---|
18 | USE traadv_muscl ! MUSCL scheme (tra_adv_muscl routine) |
---|
19 | USE traadv_muscl2 ! MUSCL2 scheme (tra_adv_muscl2 routine) |
---|
20 | USE traadv_ubs ! UBS scheme (tra_adv_ubs routine) |
---|
21 | USE traadv_qck ! QUICKEST scheme (tra_adv_qck routine) |
---|
22 | USE traadv_eiv ! eddy induced velocity (tra_adv_eiv routine) |
---|
23 | USE cla ! cross land advection (cla_traadv routine) |
---|
24 | USE ldftra_oce ! lateral diffusion coefficient on tracers |
---|
25 | USE in_out_manager ! I/O manager |
---|
26 | USE iom ! I/O module |
---|
27 | USE prtctl ! Print control |
---|
28 | USE lib_mpp ! MPP library |
---|
29 | |
---|
30 | IMPLICIT NONE |
---|
31 | PRIVATE |
---|
32 | |
---|
33 | PUBLIC tra_adv ! routine called by step module |
---|
34 | PUBLIC tra_adv_init ! routine called by opa module |
---|
35 | |
---|
36 | ! !!* Namelist namtra_adv * |
---|
37 | LOGICAL :: ln_traadv_cen2 = .TRUE. ! 2nd order centered scheme flag |
---|
38 | LOGICAL :: ln_traadv_tvd = .FALSE. ! TVD scheme flag |
---|
39 | LOGICAL :: ln_traadv_muscl = .FALSE. ! MUSCL scheme flag |
---|
40 | LOGICAL :: ln_traadv_muscl2 = .FALSE. ! MUSCL2 scheme flag |
---|
41 | LOGICAL :: ln_traadv_ubs = .FALSE. ! UBS scheme flag |
---|
42 | LOGICAL :: ln_traadv_qck = .FALSE. ! QUICKEST scheme flag |
---|
43 | |
---|
44 | INTEGER :: nadv ! choice of the type of advection scheme |
---|
45 | |
---|
46 | !! * Control permutation of array indices |
---|
47 | # include "oce_ftrans.h90" |
---|
48 | # include "dom_oce_ftrans.h90" |
---|
49 | # include "ldftra_oce_ftrans.h90" |
---|
50 | |
---|
51 | !! * Substitutions |
---|
52 | # include "domzgr_substitute.h90" |
---|
53 | # include "vectopt_loop_substitute.h90" |
---|
54 | !!---------------------------------------------------------------------- |
---|
55 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
56 | !! $Id$ |
---|
57 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
58 | !!---------------------------------------------------------------------- |
---|
59 | CONTAINS |
---|
60 | |
---|
61 | SUBROUTINE tra_adv( kt ) |
---|
62 | !!---------------------------------------------------------------------- |
---|
63 | !! *** ROUTINE tra_adv *** |
---|
64 | !! |
---|
65 | !! ** Purpose : compute the ocean tracer advection trend. |
---|
66 | !! |
---|
67 | !! ** Method : - Update (ua,va) with the advection term following nadv |
---|
68 | !!---------------------------------------------------------------------- |
---|
69 | USE wrk_nemo, ONLY: wrk_in_use, wrk_not_released |
---|
70 | USE wrk_nemo, ONLY: zun => wrk_3d_1 , zvn => wrk_3d_2 , zwn => wrk_3d_3 ! 3D workspace |
---|
71 | |
---|
72 | !! DCSE_NEMO: need additional directives for renamed module variables |
---|
73 | !FTRANS zun zvn zwn :I :I :z |
---|
74 | |
---|
75 | ! |
---|
76 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
77 | ! |
---|
78 | INTEGER :: ji, jj, jk ! dummy loop index |
---|
79 | !!---------------------------------------------------------------------- |
---|
80 | ! |
---|
81 | IF( wrk_in_use(3, 1,2,3) ) THEN |
---|
82 | CALL ctl_stop('tra_adv: requested workspace arrays unavailable') ; RETURN |
---|
83 | ENDIF |
---|
84 | ! ! set time step |
---|
85 | IF( neuler == 0 .AND. kt == nit000 ) THEN ! at nit000 |
---|
86 | r2dtra(:) = rdttra(:) ! = rdtra (restarting with Euler time stepping) |
---|
87 | ELSEIF( kt <= nit000 + 1) THEN ! at nit000 or nit000+1 |
---|
88 | r2dtra(:) = 2._wp * rdttra(:) ! = 2 rdttra (leapfrog) |
---|
89 | ENDIF |
---|
90 | ! |
---|
91 | IF( nn_cla == 1 ) CALL cla_traadv( kt ) !== Cross Land Advection ==! (hor. advection) |
---|
92 | ! |
---|
93 | ! !== effective transport ==! |
---|
94 | #if defined key_z_first |
---|
95 | DO jj = 1, jpj |
---|
96 | DO ji = 1, jpi |
---|
97 | DO jk = 1, jpkm1 |
---|
98 | zun(ji,jj,jk) = e2u(ji,jj) * fse3u(ji,jj,jk) * un(ji,jj,jk) ! eulerian transport only |
---|
99 | zvn(ji,jj,jk) = e1v(ji,jj) * fse3v(ji,jj,jk) * vn(ji,jj,jk) |
---|
100 | zwn(ji,jj,jk) = e1t(ji,jj) * e2t(ji,jj) * wn(ji,jj,jk) |
---|
101 | END DO |
---|
102 | zun(ji,jj,jpk) = 0._wp ! no transport trough the bottom |
---|
103 | zvn(ji,jj,jpk) = 0._wp ! no transport trough the bottom |
---|
104 | zwn(ji,jj,jpk) = 0._wp ! no transport trough the bottom |
---|
105 | END DO |
---|
106 | END DO |
---|
107 | #else |
---|
108 | DO jk = 1, jpkm1 |
---|
109 | zun(:,:,jk) = e2u(:,:) * fse3u(:,:,jk) * un(:,:,jk) ! eulerian transport only |
---|
110 | zvn(:,:,jk) = e1v(:,:) * fse3v(:,:,jk) * vn(:,:,jk) |
---|
111 | zwn(:,:,jk) = e1t(:,:) * e2t(:,:) * wn(:,:,jk) |
---|
112 | END DO |
---|
113 | zun(:,:,jpk) = 0._wp ! no transport trough the bottom |
---|
114 | zvn(:,:,jpk) = 0._wp ! no transport trough the bottom |
---|
115 | zwn(:,:,jpk) = 0._wp ! no transport trough the bottom |
---|
116 | #endif |
---|
117 | ! |
---|
118 | IF( lk_traldf_eiv .AND. .NOT. ln_traldf_grif ) & |
---|
119 | & CALL tra_adv_eiv( kt, zun, zvn, zwn, 'TRA' ) ! add the eiv transport (if necessary) |
---|
120 | ! |
---|
121 | CALL iom_put( "uocetr_eff", zun ) ! output effective transport |
---|
122 | CALL iom_put( "vocetr_eff", zvn ) |
---|
123 | CALL iom_put( "wocetr_eff", zwn ) |
---|
124 | |
---|
125 | SELECT CASE ( nadv ) !== compute advection trend and add it to general trend ==! |
---|
126 | CASE ( 1 ) ; CALL tra_adv_cen2 ( kt, 'TRA', zun, zvn, zwn, tsb, tsn, tsa, jpts ) ! 2nd order centered |
---|
127 | CASE ( 2 ) ; CALL tra_adv_tvd ( kt, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts ) ! TVD |
---|
128 | CASE ( 3 ) ; CALL tra_adv_muscl ( kt, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsa, jpts ) ! MUSCL |
---|
129 | CASE ( 4 ) ; CALL tra_adv_muscl2( kt, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts ) ! MUSCL2 |
---|
130 | CASE ( 5 ) ; CALL tra_adv_ubs ( kt, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts ) ! UBS |
---|
131 | CASE ( 6 ) ; CALL tra_adv_qck ( kt, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts ) ! QUICKEST |
---|
132 | ! |
---|
133 | CASE (-1 ) !== esopa: test all possibility with control print ==! |
---|
134 | CALL tra_adv_cen2 ( kt, 'TRA', zun, zvn, zwn, tsb, tsn, tsa, jpts ) |
---|
135 | CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv0 - Ta: ', mask1=tmask, & |
---|
136 | & tab3d_2=tsa(:,:,:,jp_sal), clinfo2= ' Sa: ', mask2=tmask, clinfo3='tra' ) |
---|
137 | CALL tra_adv_tvd ( kt, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts ) |
---|
138 | CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv1 - Ta: ', mask1=tmask, & |
---|
139 | & tab3d_2=tsa(:,:,:,jp_sal), clinfo2= ' Sa: ', mask2=tmask, clinfo3='tra' ) |
---|
140 | CALL tra_adv_muscl ( kt, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsa, jpts ) |
---|
141 | CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv3 - Ta: ', mask1=tmask, & |
---|
142 | & tab3d_2=tsa(:,:,:,jp_sal), clinfo2= ' Sa: ', mask2=tmask, clinfo3='tra' ) |
---|
143 | CALL tra_adv_muscl2( kt, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts ) |
---|
144 | CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv4 - Ta: ', mask1=tmask, & |
---|
145 | & tab3d_2=tsa(:,:,:,jp_sal), clinfo2= ' Sa: ', mask2=tmask, clinfo3='tra' ) |
---|
146 | CALL tra_adv_ubs ( kt, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts ) |
---|
147 | CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv5 - Ta: ', mask1=tmask, & |
---|
148 | & tab3d_2=tsa(:,:,:,jp_sal), clinfo2= ' Sa: ', mask2=tmask, clinfo3='tra' ) |
---|
149 | CALL tra_adv_qck ( kt, 'TRA', r2dtra, zun, zvn, zwn, tsb, tsn, tsa, jpts ) |
---|
150 | CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv6 - Ta: ', mask1=tmask, & |
---|
151 | & tab3d_2=tsa(:,:,:,jp_sal), clinfo2= ' Sa: ', mask2=tmask, clinfo3='tra' ) |
---|
152 | END SELECT |
---|
153 | ! |
---|
154 | ! ! print mean trends (used for debugging) |
---|
155 | IF(ln_ctl) CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' adv - Ta: ', mask1=tmask, & |
---|
156 | & tab3d_2=tsa(:,:,:,jp_sal), clinfo2= ' Sa: ', mask2=tmask, clinfo3='tra' ) |
---|
157 | ! |
---|
158 | IF( wrk_not_released(3,1,2,3) ) CALL ctl_stop('tra_adv: failed to release workspace arrays') |
---|
159 | ! |
---|
160 | END SUBROUTINE tra_adv |
---|
161 | |
---|
162 | |
---|
163 | SUBROUTINE tra_adv_init |
---|
164 | !!--------------------------------------------------------------------- |
---|
165 | !! *** ROUTINE tra_adv_init *** |
---|
166 | !! |
---|
167 | !! ** Purpose : Control the consistency between namelist options for |
---|
168 | !! tracer advection schemes and set nadv |
---|
169 | !!---------------------------------------------------------------------- |
---|
170 | INTEGER :: ioptio |
---|
171 | !! |
---|
172 | NAMELIST/namtra_adv/ ln_traadv_cen2 , ln_traadv_tvd, & |
---|
173 | & ln_traadv_muscl, ln_traadv_muscl2, & |
---|
174 | & ln_traadv_ubs , ln_traadv_qck |
---|
175 | !!---------------------------------------------------------------------- |
---|
176 | |
---|
177 | REWIND( numnam ) ! Read Namelist namtra_adv : tracer advection scheme |
---|
178 | READ ( numnam, namtra_adv ) |
---|
179 | |
---|
180 | IF(lwp) THEN ! Namelist print |
---|
181 | WRITE(numout,*) |
---|
182 | WRITE(numout,*) 'tra_adv_init : choice/control of the tracer advection scheme' |
---|
183 | WRITE(numout,*) '~~~~~~~~~~~' |
---|
184 | WRITE(numout,*) ' Namelist namtra_adv : chose a advection scheme for tracers' |
---|
185 | WRITE(numout,*) ' 2nd order advection scheme ln_traadv_cen2 = ', ln_traadv_cen2 |
---|
186 | WRITE(numout,*) ' TVD advection scheme ln_traadv_tvd = ', ln_traadv_tvd |
---|
187 | WRITE(numout,*) ' MUSCL advection scheme ln_traadv_muscl = ', ln_traadv_muscl |
---|
188 | WRITE(numout,*) ' MUSCL2 advection scheme ln_traadv_muscl2 = ', ln_traadv_muscl2 |
---|
189 | WRITE(numout,*) ' UBS advection scheme ln_traadv_ubs = ', ln_traadv_ubs |
---|
190 | WRITE(numout,*) ' QUICKEST advection scheme ln_traadv_qck = ', ln_traadv_qck |
---|
191 | ENDIF |
---|
192 | |
---|
193 | ioptio = 0 ! Parameter control |
---|
194 | IF( ln_traadv_cen2 ) ioptio = ioptio + 1 |
---|
195 | IF( ln_traadv_tvd ) ioptio = ioptio + 1 |
---|
196 | IF( ln_traadv_muscl ) ioptio = ioptio + 1 |
---|
197 | IF( ln_traadv_muscl2 ) ioptio = ioptio + 1 |
---|
198 | IF( ln_traadv_ubs ) ioptio = ioptio + 1 |
---|
199 | IF( ln_traadv_qck ) ioptio = ioptio + 1 |
---|
200 | IF( lk_esopa ) ioptio = 1 |
---|
201 | |
---|
202 | IF( ioptio /= 1 ) CALL ctl_stop( 'Choose ONE advection scheme in namelist namtra_adv' ) |
---|
203 | |
---|
204 | ! ! Set nadv |
---|
205 | IF( ln_traadv_cen2 ) nadv = 1 |
---|
206 | IF( ln_traadv_tvd ) nadv = 2 |
---|
207 | IF( ln_traadv_muscl ) nadv = 3 |
---|
208 | IF( ln_traadv_muscl2 ) nadv = 4 |
---|
209 | IF( ln_traadv_ubs ) nadv = 5 |
---|
210 | IF( ln_traadv_qck ) nadv = 6 |
---|
211 | IF( lk_esopa ) nadv = -1 |
---|
212 | |
---|
213 | IF(lwp) THEN ! Print the choice |
---|
214 | WRITE(numout,*) |
---|
215 | IF( nadv == 1 ) WRITE(numout,*) ' 2nd order scheme is used' |
---|
216 | IF( nadv == 2 ) WRITE(numout,*) ' TVD scheme is used' |
---|
217 | IF( nadv == 3 ) WRITE(numout,*) ' MUSCL scheme is used' |
---|
218 | IF( nadv == 4 ) WRITE(numout,*) ' MUSCL2 scheme is used' |
---|
219 | IF( nadv == 5 ) WRITE(numout,*) ' UBS scheme is used' |
---|
220 | IF( nadv == 6 ) WRITE(numout,*) ' QUICKEST scheme is used' |
---|
221 | IF( nadv == -1 ) WRITE(numout,*) ' esopa test: use all advection scheme' |
---|
222 | ENDIF |
---|
223 | ! |
---|
224 | END SUBROUTINE tra_adv_init |
---|
225 | |
---|
226 | !!====================================================================== |
---|
227 | END MODULE traadv |
---|