1 | MODULE istate |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE istate *** |
---|
4 | !! Ocean state : initial state setting |
---|
5 | !!===================================================================== |
---|
6 | !! History : OPA ! 1989-12 (P. Andrich) Original code |
---|
7 | !! 5.0 ! 1991-11 (G. Madec) rewritting |
---|
8 | !! 6.0 ! 1996-01 (G. Madec) terrain following coordinates |
---|
9 | !! 8.0 ! 2001-09 (M. Levy, M. Ben Jelloul) istate_eel |
---|
10 | !! 8.0 ! 2001-09 (M. Levy, M. Ben Jelloul) istate_uvg |
---|
11 | !! NEMO 1.0 ! 2003-08 (G. Madec, C. Talandier) F90: Free form, modules + EEL R5 |
---|
12 | !! - ! 2004-05 (A. Koch-Larrouy) istate_gyre |
---|
13 | !! 2.0 ! 2006-07 (S. Masson) distributed restart using iom |
---|
14 | !! 3.3 ! 2010-10 (C. Ethe) merge TRC-TRA |
---|
15 | !! 3.4 ! 2011-04 (G. Madec) Merge of dtatem and dtasal & suppression of tb,tn/sb,sn |
---|
16 | !!---------------------------------------------------------------------- |
---|
17 | |
---|
18 | !!---------------------------------------------------------------------- |
---|
19 | !! istate_init : initial state setting |
---|
20 | !! istate_tem : analytical profile for initial Temperature |
---|
21 | !! istate_sal : analytical profile for initial Salinity |
---|
22 | !! istate_eel : initial state setting of EEL R5 configuration |
---|
23 | !! istate_gyre : initial state setting of GYRE configuration |
---|
24 | !! istate_uvg : initial velocity in geostropic balance |
---|
25 | !!---------------------------------------------------------------------- |
---|
26 | USE oce ! ocean dynamics and active tracers |
---|
27 | USE dom_oce ! ocean space and time domain |
---|
28 | USE daymod ! calendar |
---|
29 | USE eosbn2 ! eq. of state, Brunt Vaisala frequency (eos routine) |
---|
30 | USE ldftra_oce ! ocean active tracers: lateral physics |
---|
31 | USE zdf_oce ! ocean vertical physics |
---|
32 | USE phycst ! physical constants |
---|
33 | USE dtatsd ! data temperature and salinity (dta_tsd routine) |
---|
34 | USE in_out_manager ! I/O manager |
---|
35 | USE iom ! I/O library |
---|
36 | USE zpshde ! partial step: hor. derivative (zps_hde routine) |
---|
37 | USE eosbn2 ! equation of state (eos bn2 routine) |
---|
38 | USE domvvl ! varying vertical mesh |
---|
39 | USE dynspg_oce ! pressure gradient schemes |
---|
40 | USE dynspg_flt ! pressure gradient schemes |
---|
41 | USE dynspg_exp ! pressure gradient schemes |
---|
42 | USE dynspg_ts ! pressure gradient schemes |
---|
43 | USE sol_oce ! ocean solver variables |
---|
44 | USE lib_mpp ! MPP library |
---|
45 | USE restart ! restart |
---|
46 | USE wrk_nemo ! Memory allocation |
---|
47 | USE timing ! Timing |
---|
48 | |
---|
49 | IMPLICIT NONE |
---|
50 | PRIVATE |
---|
51 | |
---|
52 | PUBLIC istate_init ! routine called by step.F90 |
---|
53 | |
---|
54 | !! * Substitutions |
---|
55 | # include "domzgr_substitute.h90" |
---|
56 | # include "vectopt_loop_substitute.h90" |
---|
57 | !!---------------------------------------------------------------------- |
---|
58 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
59 | !! $Id$ |
---|
60 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
61 | !!---------------------------------------------------------------------- |
---|
62 | CONTAINS |
---|
63 | |
---|
64 | SUBROUTINE istate_init |
---|
65 | !!---------------------------------------------------------------------- |
---|
66 | !! *** ROUTINE istate_init *** |
---|
67 | !! |
---|
68 | !! ** Purpose : Initialization of the dynamics and tracer fields. |
---|
69 | !!---------------------------------------------------------------------- |
---|
70 | ! - ML - needed for initialization of e3t_b |
---|
71 | INTEGER :: jk ! dummy loop indice |
---|
72 | !!---------------------------------------------------------------------- |
---|
73 | ! |
---|
74 | IF( nn_timing == 1 ) CALL timing_start('istate_init') |
---|
75 | ! |
---|
76 | |
---|
77 | IF(lwp) WRITE(numout,*) |
---|
78 | IF(lwp) WRITE(numout,*) 'istate_ini : Initialization of the dynamics and tracers' |
---|
79 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~' |
---|
80 | |
---|
81 | CALL dta_tsd_init ! Initialisation of T & S input data |
---|
82 | |
---|
83 | rhd (:,:,: ) = 0.e0 |
---|
84 | rhop (:,:,: ) = 0.e0 |
---|
85 | rn2 (:,:,: ) = 0.e0 |
---|
86 | tsa (:,:,:,:) = 0.e0 |
---|
87 | |
---|
88 | IF( ln_rstart ) THEN ! Restart from a file |
---|
89 | ! ! ------------------- |
---|
90 | neuler = 1 ! Set time-step indicator at nit000 (leap-frog) |
---|
91 | CALL rst_read ! Read the restart file |
---|
92 | ! ! define e3u_b, e3v_b from e3t_b read in restart file |
---|
93 | CALL dom_vvl_2( nit000, fse3u_b(:,:,:), fse3v_b(:,:,:) ) |
---|
94 | CALL day_init ! model calendar (using both namelist and restart infos) |
---|
95 | ELSE |
---|
96 | ! ! Start from rest |
---|
97 | ! ! --------------- |
---|
98 | numror = 0 ! define numror = 0 -> no restart file to read |
---|
99 | neuler = 0 ! Set time-step indicator at nit000 (euler forward) |
---|
100 | CALL day_init ! model calendar (using both namelist and restart infos) |
---|
101 | ! ! Initialization of ocean to zero |
---|
102 | ! before fields ! now fields |
---|
103 | sshb (:,:) = 0._wp ; sshn (:,:) = 0._wp |
---|
104 | ub (:,:,:) = 0._wp ; un (:,:,:) = 0._wp |
---|
105 | vb (:,:,:) = 0._wp ; vn (:,:,:) = 0._wp |
---|
106 | rotb (:,:,:) = 0._wp ; rotn (:,:,:) = 0._wp |
---|
107 | hdivb(:,:,:) = 0._wp ; hdivn(:,:,:) = 0._wp |
---|
108 | ! |
---|
109 | IF( cp_cfg == 'eel' ) THEN |
---|
110 | CALL istate_eel ! EEL configuration : start from pre-defined U,V T-S fields |
---|
111 | ELSEIF( cp_cfg == 'gyre' ) THEN |
---|
112 | CALL istate_gyre ! GYRE configuration : start from pre-defined T-S fields |
---|
113 | ELSEIF( ln_tsd_init ) THEN ! Initial T-S fields read in files |
---|
114 | CALL dta_tsd( nit000, tsb ) ! read 3D T and S data at nit000 |
---|
115 | tsn(:,:,:,:) = tsb(:,:,:,:) |
---|
116 | ! |
---|
117 | ELSE ! Initial T-S fields defined analytically |
---|
118 | CALL istate_t_s |
---|
119 | ENDIF |
---|
120 | ! |
---|
121 | CALL eos( tsb, rhd, rhop ) ! before potential and in situ densities |
---|
122 | #if ! defined key_c1d |
---|
123 | IF( ln_zps ) CALL zps_hde( nit000, jpts, tsb, gtsu, gtsv, & ! zps: before hor. gradient |
---|
124 | & rhd, gru , grv ) ! of t,s,rd at ocean bottom |
---|
125 | #endif |
---|
126 | ! |
---|
127 | ! - ML - sshn could be modified by istate_eel, so that initialization of fse3t_b is done here |
---|
128 | IF( lk_vvl ) THEN |
---|
129 | DO jk = 1, jpk |
---|
130 | fse3t_b(:,:,jk) = fse3t_n(:,:,jk) |
---|
131 | ENDDO |
---|
132 | ENDIF |
---|
133 | ! ! define e3u_b, e3v_b from e3t_b initialized in domzgr |
---|
134 | CALL dom_vvl_2( nit000, fse3u_b(:,:,:), fse3v_b(:,:,:) ) |
---|
135 | ! |
---|
136 | ENDIF |
---|
137 | ! |
---|
138 | IF( lk_agrif ) THEN ! read free surface arrays in restart file |
---|
139 | IF( ln_rstart ) THEN |
---|
140 | IF( lk_dynspg_flt ) THEN ! read or initialize the following fields |
---|
141 | ! ! gcx, gcxb for agrif_opa_init |
---|
142 | IF( sol_oce_alloc() > 0 ) CALL ctl_stop('agrif sol_oce_alloc: allocation of arrays failed') |
---|
143 | CALL flt_rst( nit000, 'READ' ) |
---|
144 | ENDIF |
---|
145 | ENDIF ! explicit case not coded yet with AGRIF |
---|
146 | ENDIF |
---|
147 | ! |
---|
148 | IF( nn_timing == 1 ) CALL timing_stop('istate_init') |
---|
149 | ! |
---|
150 | END SUBROUTINE istate_init |
---|
151 | |
---|
152 | SUBROUTINE istate_t_s |
---|
153 | !!--------------------------------------------------------------------- |
---|
154 | !! *** ROUTINE istate_t_s *** |
---|
155 | !! |
---|
156 | !! ** Purpose : Intialization of the temperature field with an |
---|
157 | !! analytical profile or a file (i.e. in EEL configuration) |
---|
158 | !! |
---|
159 | !! ** Method : - temperature: use Philander analytic profile |
---|
160 | !! - salinity : use to a constant value 35.5 |
---|
161 | !! |
---|
162 | !! References : Philander ??? |
---|
163 | !!---------------------------------------------------------------------- |
---|
164 | INTEGER :: ji, jj, jk |
---|
165 | REAL(wp) :: zsal = 35.50 |
---|
166 | !!---------------------------------------------------------------------- |
---|
167 | ! |
---|
168 | IF(lwp) WRITE(numout,*) |
---|
169 | IF(lwp) WRITE(numout,*) 'istate_t_s : Philander s initial temperature profile' |
---|
170 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~ and constant salinity (',zsal,' psu)' |
---|
171 | ! |
---|
172 | DO jk = 1, jpk |
---|
173 | tsn(:,:,jk,jp_tem) = ( ( ( 7.5 - 0. * ABS( gphit(:,:) )/30. ) * ( 1.-TANH((fsdept(:,:,jk)-80.)/30.) ) & |
---|
174 | & + 10. * ( 5000. - fsdept(:,:,jk) ) /5000.) ) * tmask(:,:,jk) |
---|
175 | tsb(:,:,jk,jp_tem) = tsn(:,:,jk,jp_tem) |
---|
176 | END DO |
---|
177 | tsn(:,:,:,jp_sal) = zsal * tmask(:,:,:) |
---|
178 | tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal) |
---|
179 | ! |
---|
180 | END SUBROUTINE istate_t_s |
---|
181 | |
---|
182 | SUBROUTINE istate_eel |
---|
183 | !!---------------------------------------------------------------------- |
---|
184 | !! *** ROUTINE istate_eel *** |
---|
185 | !! |
---|
186 | !! ** Purpose : Initialization of the dynamics and tracers for EEL R5 |
---|
187 | !! configuration (channel with or without a topographic bump) |
---|
188 | !! |
---|
189 | !! ** Method : - set temprature field |
---|
190 | !! - set salinity field |
---|
191 | !! - set velocity field including horizontal divergence |
---|
192 | !! and relative vorticity fields |
---|
193 | !!---------------------------------------------------------------------- |
---|
194 | USE divcur ! hor. divergence & rel. vorticity (div_cur routine) |
---|
195 | USE iom |
---|
196 | |
---|
197 | INTEGER :: inum ! temporary logical unit |
---|
198 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
199 | INTEGER :: ijloc |
---|
200 | REAL(wp) :: zh1, zh2, zslope, zcst, zfcor ! temporary scalars |
---|
201 | REAL(wp) :: zt1 = 15._wp ! surface temperature value (EEL R5) |
---|
202 | REAL(wp) :: zt2 = 5._wp ! bottom temperature value (EEL R5) |
---|
203 | REAL(wp) :: zsal = 35.0_wp ! constant salinity (EEL R2, R5 and R6) |
---|
204 | REAL(wp) :: zueel = 0.1_wp ! constant uniform zonal velocity (EEL R5) |
---|
205 | REAL(wp), DIMENSION(jpiglo,jpjglo) :: zssh ! initial ssh over the global domain |
---|
206 | !!---------------------------------------------------------------------- |
---|
207 | |
---|
208 | SELECT CASE ( jp_cfg ) |
---|
209 | ! ! ==================== |
---|
210 | CASE ( 5 ) ! EEL R5 configuration |
---|
211 | ! ! ==================== |
---|
212 | ! |
---|
213 | ! set temperature field with a linear profile |
---|
214 | ! ------------------------------------------- |
---|
215 | IF(lwp) WRITE(numout,*) |
---|
216 | IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: linear temperature profile' |
---|
217 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~' |
---|
218 | ! |
---|
219 | zh1 = gdept_0( 1 ) |
---|
220 | zh2 = gdept_0(jpkm1) |
---|
221 | ! |
---|
222 | zslope = ( zt1 - zt2 ) / ( zh1 - zh2 ) |
---|
223 | zcst = ( zt1 * ( zh1 - zh2) - ( zt1 - zt2 ) * zh1 ) / ( zh1 - zh2 ) |
---|
224 | ! |
---|
225 | DO jk = 1, jpk |
---|
226 | tsn(:,:,jk,jp_tem) = ( zt2 + zt1 * exp( - fsdept(:,:,jk) / 1000 ) ) * tmask(:,:,jk) |
---|
227 | tsb(:,:,jk,jp_tem) = tsn(:,:,jk,jp_tem) |
---|
228 | END DO |
---|
229 | ! |
---|
230 | IF(lwp) CALL prizre( tsn(:,:,:,jp_tem), jpi , jpj , jpk , jpj/2 , & |
---|
231 | & 1 , jpi , 5 , 1 , jpk , & |
---|
232 | & 1 , 1. , numout ) |
---|
233 | ! |
---|
234 | ! set salinity field to a constant value |
---|
235 | ! -------------------------------------- |
---|
236 | IF(lwp) WRITE(numout,*) |
---|
237 | IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal |
---|
238 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~' |
---|
239 | ! |
---|
240 | tsn(:,:,:,jp_sal) = zsal * tmask(:,:,:) |
---|
241 | tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal) |
---|
242 | ! |
---|
243 | ! set the dynamics: U,V, hdiv, rot (and ssh if necessary) |
---|
244 | ! ---------------- |
---|
245 | ! Start EEL5 configuration with barotropic geostrophic velocities |
---|
246 | ! according the sshb and sshn SSH imposed. |
---|
247 | ! we assume a uniform grid (hence the use of e1t(1,1) for delta_y) |
---|
248 | ! we use the Coriolis frequency at mid-channel. |
---|
249 | ub(:,:,:) = zueel * umask(:,:,:) |
---|
250 | un(:,:,:) = ub(:,:,:) |
---|
251 | ijloc = mj0(INT(jpjglo-1)/2) |
---|
252 | zfcor = ff(1,ijloc) |
---|
253 | ! |
---|
254 | DO jj = 1, jpjglo |
---|
255 | zssh(:,jj) = - (FLOAT(jj)- FLOAT(jpjglo-1)/2.)*zueel*e1t(1,1)*zfcor/grav |
---|
256 | END DO |
---|
257 | ! |
---|
258 | IF(lwp) THEN |
---|
259 | WRITE(numout,*) ' Uniform zonal velocity for EEL R5:',zueel |
---|
260 | WRITE(numout,*) ' Geostrophic SSH profile as a function of y:' |
---|
261 | WRITE(numout,'(12(1x,f6.2))') zssh(1,:) |
---|
262 | ENDIF |
---|
263 | ! |
---|
264 | DO jj = 1, nlcj |
---|
265 | DO ji = 1, nlci |
---|
266 | sshb(ji,jj) = zssh( mig(ji) , mjg(jj) ) * tmask(ji,jj,1) |
---|
267 | END DO |
---|
268 | END DO |
---|
269 | sshb(nlci+1:jpi, : ) = 0.e0 ! set to zero extra mpp columns |
---|
270 | sshb( : ,nlcj+1:jpj) = 0.e0 ! set to zero extra mpp rows |
---|
271 | ! |
---|
272 | sshn(:,:) = sshb(:,:) ! set now ssh to the before value |
---|
273 | ! |
---|
274 | IF( nn_rstssh /= 0 ) THEN |
---|
275 | nn_rstssh = 0 ! hand-made initilization of ssh |
---|
276 | CALL ctl_warn( 'istate_eel: force nn_rstssh = 0' ) |
---|
277 | ENDIF |
---|
278 | ! |
---|
279 | CALL div_cur( nit000 ) ! horizontal divergence and relative vorticity (curl) |
---|
280 | ! N.B. the vertical velocity will be computed from the horizontal divergence field |
---|
281 | ! in istate by a call to wzv routine |
---|
282 | |
---|
283 | |
---|
284 | ! ! ========================== |
---|
285 | CASE ( 2 , 6 ) ! EEL R2 or R6 configuration |
---|
286 | ! ! ========================== |
---|
287 | ! |
---|
288 | ! set temperature field with a NetCDF file |
---|
289 | ! ---------------------------------------- |
---|
290 | IF(lwp) WRITE(numout,*) |
---|
291 | IF(lwp) WRITE(numout,*) 'istate_eel : EEL R2 or R6: read initial temperature in a NetCDF file' |
---|
292 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~' |
---|
293 | ! |
---|
294 | CALL iom_open ( 'eel.initemp', inum ) |
---|
295 | CALL iom_get ( inum, jpdom_data, 'initemp', tsb(:,:,:,jp_tem) ) ! read before temprature (tb) |
---|
296 | CALL iom_close( inum ) |
---|
297 | ! |
---|
298 | tsn(:,:,:,jp_tem) = tsb(:,:,:,jp_tem) ! set nox temperature to tb |
---|
299 | ! |
---|
300 | IF(lwp) CALL prizre( tsn(:,:,:,jp_tem), jpi , jpj , jpk , jpj/2 , & |
---|
301 | & 1 , jpi , 5 , 1 , jpk , & |
---|
302 | & 1 , 1. , numout ) |
---|
303 | ! |
---|
304 | ! set salinity field to a constant value |
---|
305 | ! -------------------------------------- |
---|
306 | IF(lwp) WRITE(numout,*) |
---|
307 | IF(lwp) WRITE(numout,*) 'istate_eel : EEL R5: constant salinity field, S = ', zsal |
---|
308 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~' |
---|
309 | ! |
---|
310 | tsn(:,:,:,jp_sal) = zsal * tmask(:,:,:) |
---|
311 | tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal) |
---|
312 | ! |
---|
313 | ! ! =========================== |
---|
314 | CASE DEFAULT ! NONE existing configuration |
---|
315 | ! ! =========================== |
---|
316 | WRITE(ctmp1,*) 'EEL with a ', jp_cfg,' km resolution is not coded' |
---|
317 | CALL ctl_stop( ctmp1 ) |
---|
318 | ! |
---|
319 | END SELECT |
---|
320 | ! |
---|
321 | END SUBROUTINE istate_eel |
---|
322 | |
---|
323 | |
---|
324 | SUBROUTINE istate_gyre |
---|
325 | !!---------------------------------------------------------------------- |
---|
326 | !! *** ROUTINE istate_gyre *** |
---|
327 | !! |
---|
328 | !! ** Purpose : Initialization of the dynamics and tracers for GYRE |
---|
329 | !! configuration (double gyre with rotated domain) |
---|
330 | !! |
---|
331 | !! ** Method : - set temprature field |
---|
332 | !! - set salinity field |
---|
333 | !!---------------------------------------------------------------------- |
---|
334 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
335 | INTEGER :: inum ! temporary logical unit |
---|
336 | INTEGER, PARAMETER :: ntsinit = 0 ! (0/1) (analytical/input data files) T&S initialization |
---|
337 | !!---------------------------------------------------------------------- |
---|
338 | |
---|
339 | SELECT CASE ( ntsinit) |
---|
340 | |
---|
341 | CASE ( 0 ) ! analytical T/S profil deduced from LEVITUS |
---|
342 | IF(lwp) WRITE(numout,*) |
---|
343 | IF(lwp) WRITE(numout,*) 'istate_gyre : initial analytical T and S profil deduced from LEVITUS ' |
---|
344 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~~' |
---|
345 | |
---|
346 | DO jk = 1, jpk |
---|
347 | DO jj = 1, jpj |
---|
348 | DO ji = 1, jpi |
---|
349 | tsn(ji,jj,jk,jp_tem) = ( 16. - 12. * TANH( (fsdept(ji,jj,jk) - 400) / 700 ) ) & |
---|
350 | & * (-TANH( (500-fsdept(ji,jj,jk)) / 150 ) + 1) / 2 & |
---|
351 | & + ( 15. * ( 1. - TANH( (fsdept(ji,jj,jk)-50.) / 1500.) ) & |
---|
352 | & - 1.4 * TANH((fsdept(ji,jj,jk)-100.) / 100.) & |
---|
353 | & + 7. * (1500. - fsdept(ji,jj,jk)) / 1500. ) & |
---|
354 | & * (-TANH( (fsdept(ji,jj,jk) - 500) / 150) + 1) / 2 |
---|
355 | tsn(ji,jj,jk,jp_tem) = tsn(ji,jj,jk,jp_tem) * tmask(ji,jj,jk) |
---|
356 | tsb(ji,jj,jk,jp_tem) = tsn(ji,jj,jk,jp_tem) |
---|
357 | |
---|
358 | tsn(ji,jj,jk,jp_sal) = ( 36.25 - 1.13 * TANH( (fsdept(ji,jj,jk) - 305) / 460 ) ) & |
---|
359 | & * (-TANH((500 - fsdept(ji,jj,jk)) / 150) + 1) / 2 & |
---|
360 | & + ( 35.55 + 1.25 * (5000. - fsdept(ji,jj,jk)) / 5000. & |
---|
361 | & - 1.62 * TANH( (fsdept(ji,jj,jk) - 60. ) / 650. ) & |
---|
362 | & + 0.2 * TANH( (fsdept(ji,jj,jk) - 35. ) / 100. ) & |
---|
363 | & + 0.2 * TANH( (fsdept(ji,jj,jk) - 1000.) / 5000.) ) & |
---|
364 | & * (-TANH((fsdept(ji,jj,jk) - 500) / 150) + 1) / 2 |
---|
365 | tsn(ji,jj,jk,jp_sal) = tsn(ji,jj,jk,jp_sal) * tmask(ji,jj,jk) |
---|
366 | tsb(ji,jj,jk,jp_sal) = tsn(ji,jj,jk,jp_sal) |
---|
367 | END DO |
---|
368 | END DO |
---|
369 | END DO |
---|
370 | |
---|
371 | CASE ( 1 ) ! T/S data fields read in dta_tem.nc/data_sal.nc files |
---|
372 | IF(lwp) WRITE(numout,*) |
---|
373 | IF(lwp) WRITE(numout,*) 'istate_gyre : initial T and S read from dta_tem.nc/data_sal.nc files' |
---|
374 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~~' |
---|
375 | IF(lwp) WRITE(numout,*) ' NetCDF FORMAT' |
---|
376 | |
---|
377 | ! Read temperature field |
---|
378 | ! ---------------------- |
---|
379 | CALL iom_open ( 'data_tem', inum ) |
---|
380 | CALL iom_get ( inum, jpdom_data, 'votemper', tsn(:,:,:,jp_tem) ) |
---|
381 | CALL iom_close( inum ) |
---|
382 | |
---|
383 | tsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem) * tmask(:,:,:) |
---|
384 | tsb(:,:,:,jp_tem) = tsn(:,:,:,jp_tem) |
---|
385 | |
---|
386 | ! Read salinity field |
---|
387 | ! ------------------- |
---|
388 | CALL iom_open ( 'data_sal', inum ) |
---|
389 | CALL iom_get ( inum, jpdom_data, 'vosaline', tsn(:,:,:,jp_sal) ) |
---|
390 | CALL iom_close( inum ) |
---|
391 | |
---|
392 | tsn(:,:,:,jp_sal) = tsn(:,:,:,jp_sal) * tmask(:,:,:) |
---|
393 | tsb(:,:,:,jp_sal) = tsn(:,:,:,jp_sal) |
---|
394 | |
---|
395 | END SELECT |
---|
396 | |
---|
397 | IF(lwp) THEN |
---|
398 | WRITE(numout,*) |
---|
399 | WRITE(numout,*) ' Initial temperature and salinity profiles:' |
---|
400 | WRITE(numout, "(9x,' level gdept_0 temperature salinity ')" ) |
---|
401 | WRITE(numout, "(10x, i4, 3f10.2)" ) ( jk, gdept_0(jk), tsn(2,2,jk,jp_tem), tsn(2,2,jk,jp_sal), jk = 1, jpk ) |
---|
402 | ENDIF |
---|
403 | |
---|
404 | END SUBROUTINE istate_gyre |
---|
405 | |
---|
406 | SUBROUTINE istate_uvg |
---|
407 | !!---------------------------------------------------------------------- |
---|
408 | !! *** ROUTINE istate_uvg *** |
---|
409 | !! |
---|
410 | !! ** Purpose : Compute the geostrophic velocities from (tn,sn) fields |
---|
411 | !! |
---|
412 | !! ** Method : Using the hydrostatic hypothesis the now hydrostatic |
---|
413 | !! pressure is computed by integrating the in-situ density from the |
---|
414 | !! surface to the bottom. |
---|
415 | !! p=integral [ rau*g dz ] |
---|
416 | !!---------------------------------------------------------------------- |
---|
417 | USE dynspg ! surface pressure gradient (dyn_spg routine) |
---|
418 | USE divcur ! hor. divergence & rel. vorticity (div_cur routine) |
---|
419 | USE lbclnk ! ocean lateral boundary condition (or mpp link) |
---|
420 | |
---|
421 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
422 | INTEGER :: indic ! ??? |
---|
423 | REAL(wp) :: zmsv, zphv, zmsu, zphu, zalfg ! temporary scalars |
---|
424 | REAL(wp), POINTER, DIMENSION(:,:,:) :: zprn |
---|
425 | !!---------------------------------------------------------------------- |
---|
426 | ! |
---|
427 | CALL wrk_alloc( jpi, jpj, jpk, zprn) |
---|
428 | ! |
---|
429 | IF(lwp) WRITE(numout,*) |
---|
430 | IF(lwp) WRITE(numout,*) 'istate_uvg : Start from Geostrophy' |
---|
431 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~' |
---|
432 | |
---|
433 | ! Compute the now hydrostatic pressure |
---|
434 | ! ------------------------------------ |
---|
435 | |
---|
436 | zalfg = 0.5 * grav * rau0 |
---|
437 | |
---|
438 | zprn(:,:,1) = zalfg * fse3w(:,:,1) * ( 1 + rhd(:,:,1) ) ! Surface value |
---|
439 | |
---|
440 | DO jk = 2, jpkm1 ! Vertical integration from the surface |
---|
441 | zprn(:,:,jk) = zprn(:,:,jk-1) & |
---|
442 | & + zalfg * fse3w(:,:,jk) * ( 2. + rhd(:,:,jk) + rhd(:,:,jk-1) ) |
---|
443 | END DO |
---|
444 | |
---|
445 | ! Compute geostrophic balance |
---|
446 | ! --------------------------- |
---|
447 | DO jk = 1, jpkm1 |
---|
448 | DO jj = 2, jpjm1 |
---|
449 | DO ji = fs_2, fs_jpim1 ! vertor opt. |
---|
450 | zmsv = 1. / MAX( umask(ji-1,jj+1,jk) + umask(ji ,jj+1,jk) & |
---|
451 | + umask(ji-1,jj ,jk) + umask(ji ,jj ,jk) , 1. ) |
---|
452 | zphv = ( zprn(ji ,jj+1,jk) - zprn(ji-1,jj+1,jk) ) * umask(ji-1,jj+1,jk) / e1u(ji-1,jj+1) & |
---|
453 | + ( zprn(ji+1,jj+1,jk) - zprn(ji ,jj+1,jk) ) * umask(ji ,jj+1,jk) / e1u(ji ,jj+1) & |
---|
454 | + ( zprn(ji ,jj ,jk) - zprn(ji-1,jj ,jk) ) * umask(ji-1,jj ,jk) / e1u(ji-1,jj ) & |
---|
455 | + ( zprn(ji+1,jj ,jk) - zprn(ji ,jj ,jk) ) * umask(ji ,jj ,jk) / e1u(ji ,jj ) |
---|
456 | zphv = 1. / rau0 * zphv * zmsv * vmask(ji,jj,jk) |
---|
457 | |
---|
458 | zmsu = 1. / MAX( vmask(ji+1,jj ,jk) + vmask(ji ,jj ,jk) & |
---|
459 | + vmask(ji+1,jj-1,jk) + vmask(ji ,jj-1,jk) , 1. ) |
---|
460 | zphu = ( zprn(ji+1,jj+1,jk) - zprn(ji+1,jj ,jk) ) * vmask(ji+1,jj ,jk) / e2v(ji+1,jj ) & |
---|
461 | + ( zprn(ji ,jj+1,jk) - zprn(ji ,jj ,jk) ) * vmask(ji ,jj ,jk) / e2v(ji ,jj ) & |
---|
462 | + ( zprn(ji+1,jj ,jk) - zprn(ji+1,jj-1,jk) ) * vmask(ji+1,jj-1,jk) / e2v(ji+1,jj-1) & |
---|
463 | + ( zprn(ji ,jj ,jk) - zprn(ji ,jj-1,jk) ) * vmask(ji ,jj-1,jk) / e2v(ji ,jj-1) |
---|
464 | zphu = 1. / rau0 * zphu * zmsu * umask(ji,jj,jk) |
---|
465 | |
---|
466 | ! Compute the geostrophic velocities |
---|
467 | un(ji,jj,jk) = -2. * zphu / ( ff(ji,jj) + ff(ji ,jj-1) ) |
---|
468 | vn(ji,jj,jk) = 2. * zphv / ( ff(ji,jj) + ff(ji-1,jj ) ) |
---|
469 | END DO |
---|
470 | END DO |
---|
471 | END DO |
---|
472 | |
---|
473 | IF(lwp) WRITE(numout,*) ' we force to zero bottom velocity' |
---|
474 | |
---|
475 | ! Susbtract the bottom velocity (level jpk-1 for flat bottom case) |
---|
476 | ! to have a zero bottom velocity |
---|
477 | |
---|
478 | DO jk = 1, jpkm1 |
---|
479 | un(:,:,jk) = ( un(:,:,jk) - un(:,:,jpkm1) ) * umask(:,:,jk) |
---|
480 | vn(:,:,jk) = ( vn(:,:,jk) - vn(:,:,jpkm1) ) * vmask(:,:,jk) |
---|
481 | END DO |
---|
482 | |
---|
483 | CALL lbc_lnk( un, 'U', -1. ) |
---|
484 | CALL lbc_lnk( vn, 'V', -1. ) |
---|
485 | |
---|
486 | ub(:,:,:) = un(:,:,:) |
---|
487 | vb(:,:,:) = vn(:,:,:) |
---|
488 | |
---|
489 | ! WARNING !!!!! |
---|
490 | ! after initializing u and v, we need to calculate the initial streamfunction bsf. |
---|
491 | ! Otherwise, only the trend will be computed and the model will blow up (inconsistency). |
---|
492 | ! to do that, we call dyn_spg with a special trick: |
---|
493 | ! we fill ua and va with the velocities divided by dt, and the streamfunction will be brought to the |
---|
494 | ! right value assuming the velocities have been set up in one time step. |
---|
495 | ! we then set bsfd to zero (first guess for next step is d(psi)/dt = 0.) |
---|
496 | ! sets up s false trend to calculate the barotropic streamfunction. |
---|
497 | |
---|
498 | ua(:,:,:) = ub(:,:,:) / rdt |
---|
499 | va(:,:,:) = vb(:,:,:) / rdt |
---|
500 | |
---|
501 | ! calls dyn_spg. we assume euler time step, starting from rest. |
---|
502 | indic = 0 |
---|
503 | CALL dyn_spg( nit000, indic ) ! surface pressure gradient |
---|
504 | |
---|
505 | ! the new velocity is ua*rdt |
---|
506 | |
---|
507 | CALL lbc_lnk( ua, 'U', -1. ) |
---|
508 | CALL lbc_lnk( va, 'V', -1. ) |
---|
509 | |
---|
510 | ub(:,:,:) = ua(:,:,:) * rdt |
---|
511 | vb(:,:,:) = va(:,:,:) * rdt |
---|
512 | ua(:,:,:) = 0.e0 |
---|
513 | va(:,:,:) = 0.e0 |
---|
514 | un(:,:,:) = ub(:,:,:) |
---|
515 | vn(:,:,:) = vb(:,:,:) |
---|
516 | |
---|
517 | ! Compute the divergence and curl |
---|
518 | |
---|
519 | CALL div_cur( nit000 ) ! now horizontal divergence and curl |
---|
520 | |
---|
521 | hdivb(:,:,:) = hdivn(:,:,:) ! set the before to the now value |
---|
522 | rotb (:,:,:) = rotn (:,:,:) ! set the before to the now value |
---|
523 | ! |
---|
524 | CALL wrk_dealloc( jpi, jpj, jpk, zprn) |
---|
525 | ! |
---|
526 | END SUBROUTINE istate_uvg |
---|
527 | |
---|
528 | !!===================================================================== |
---|
529 | END MODULE istate |
---|