1 | MODULE diawri |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE diawri *** |
---|
4 | !! Ocean diagnostics : write ocean output files |
---|
5 | !!===================================================================== |
---|
6 | !! History : OPA ! 1991-03 (M.-A. Foujols) Original code |
---|
7 | !! 4.0 ! 1991-11 (G. Madec) |
---|
8 | !! ! 1992-06 (M. Imbard) correction restart file |
---|
9 | !! ! 1992-07 (M. Imbard) split into diawri and rstwri |
---|
10 | !! ! 1993-03 (M. Imbard) suppress writibm |
---|
11 | !! ! 1998-01 (C. Levy) NETCDF format using ioipsl INTERFACE |
---|
12 | !! ! 1999-02 (E. Guilyardi) name of netCDF files + variables |
---|
13 | !! 8.2 ! 2000-06 (M. Imbard) Original code (diabort.F) |
---|
14 | !! NEMO 1.0 ! 2002-06 (A.Bozec, E. Durand) Original code (diainit.F) |
---|
15 | !! - ! 2002-09 (G. Madec) F90: Free form and module |
---|
16 | !! - ! 2002-12 (G. Madec) merge of diabort and diainit, F90 |
---|
17 | !! ! 2005-11 (V. Garnier) Surface pressure gradient organization |
---|
18 | !! 3.2 ! 2008-11 (B. Lemaire) creation from old diawri |
---|
19 | !! 3.7 ! 2014-01 (G. Madec) remove eddy induced velocity from no-IOM output |
---|
20 | !! ! change name of output variables in dia_wri_state |
---|
21 | !!---------------------------------------------------------------------- |
---|
22 | |
---|
23 | !!---------------------------------------------------------------------- |
---|
24 | !! dia_wri : create the standart output files |
---|
25 | !! dia_wri_state : create an output NetCDF file for a single instantaeous ocean state and forcing fields |
---|
26 | !!---------------------------------------------------------------------- |
---|
27 | USE oce ! ocean dynamics and tracers |
---|
28 | USE dom_oce ! ocean space and time domain |
---|
29 | USE phycst ! physical constants |
---|
30 | USE dianam ! build name of file (routine) |
---|
31 | USE diahth ! thermocline diagnostics |
---|
32 | USE dynadv , ONLY: ln_dynadv_vec |
---|
33 | USE icb_oce ! Icebergs |
---|
34 | USE icbdia ! Iceberg budgets |
---|
35 | USE ldftra ! lateral physics: eddy diffusivity coef. |
---|
36 | USE ldfdyn ! lateral physics: eddy viscosity coef. |
---|
37 | USE sbc_oce ! Surface boundary condition: ocean fields |
---|
38 | USE sbc_ice ! Surface boundary condition: ice fields |
---|
39 | USE sbcssr ! restoring term toward SST/SSS climatology |
---|
40 | USE sbcwave ! wave parameters |
---|
41 | USE wet_dry ! wetting and drying |
---|
42 | USE zdf_oce ! ocean vertical physics |
---|
43 | USE zdfdrg ! ocean vertical physics: top/bottom friction |
---|
44 | USE zdfmxl ! mixed layer |
---|
45 | ! |
---|
46 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
47 | USE in_out_manager ! I/O manager |
---|
48 | USE diatmb ! Top,middle,bottom output |
---|
49 | USE dia25h ! 25h Mean output |
---|
50 | USE iom ! |
---|
51 | USE ioipsl ! |
---|
52 | |
---|
53 | #if defined key_si3 |
---|
54 | USE ice |
---|
55 | USE icewri |
---|
56 | #endif |
---|
57 | USE lib_mpp ! MPP library |
---|
58 | USE timing ! preformance summary |
---|
59 | USE diurnal_bulk ! diurnal warm layer |
---|
60 | USE cool_skin ! Cool skin |
---|
61 | |
---|
62 | IMPLICIT NONE |
---|
63 | PRIVATE |
---|
64 | |
---|
65 | PUBLIC dia_wri ! routines called by step.F90 |
---|
66 | PUBLIC dia_wri_state |
---|
67 | PUBLIC dia_wri_alloc ! Called by nemogcm module |
---|
68 | |
---|
69 | INTEGER :: nid_T, nz_T, nh_T, ndim_T, ndim_hT ! grid_T file |
---|
70 | INTEGER :: nb_T , ndim_bT ! grid_T file |
---|
71 | INTEGER :: nid_U, nz_U, nh_U, ndim_U, ndim_hU ! grid_U file |
---|
72 | INTEGER :: nid_V, nz_V, nh_V, ndim_V, ndim_hV ! grid_V file |
---|
73 | INTEGER :: nid_W, nz_W, nh_W ! grid_W file |
---|
74 | INTEGER :: ndex(1) ! ??? |
---|
75 | INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_hT, ndex_hU, ndex_hV |
---|
76 | INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_T, ndex_U, ndex_V |
---|
77 | INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_bT |
---|
78 | |
---|
79 | !! * Substitutions |
---|
80 | # include "vectopt_loop_substitute.h90" |
---|
81 | !!---------------------------------------------------------------------- |
---|
82 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
83 | !! $Id: diawri.F90 9598 2018-05-15 22:47:16Z nicolasmartin $ |
---|
84 | !! Software governed by the CeCILL licence (./LICENSE) |
---|
85 | !!---------------------------------------------------------------------- |
---|
86 | CONTAINS |
---|
87 | |
---|
88 | #if defined key_iomput |
---|
89 | !!---------------------------------------------------------------------- |
---|
90 | !! 'key_iomput' use IOM library |
---|
91 | !!---------------------------------------------------------------------- |
---|
92 | INTEGER FUNCTION dia_wri_alloc() |
---|
93 | ! |
---|
94 | dia_wri_alloc = 0 |
---|
95 | ! |
---|
96 | END FUNCTION dia_wri_alloc |
---|
97 | |
---|
98 | |
---|
99 | SUBROUTINE dia_wri( kt ) |
---|
100 | !!--------------------------------------------------------------------- |
---|
101 | !! *** ROUTINE dia_wri *** |
---|
102 | !! |
---|
103 | !! ** Purpose : Standard output of opa: dynamics and tracer fields |
---|
104 | !! NETCDF format is used by default |
---|
105 | !! |
---|
106 | !! ** Method : use iom_put |
---|
107 | !!---------------------------------------------------------------------- |
---|
108 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
109 | !! |
---|
110 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
111 | INTEGER :: ikbot ! local integer |
---|
112 | REAL(wp):: zztmp , zztmpx ! local scalar |
---|
113 | REAL(wp):: zztmp2, zztmpy ! - - |
---|
114 | REAL(wp), DIMENSION(jpi,jpj) :: z2d ! 2D workspace |
---|
115 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: z3d ! 3D workspace |
---|
116 | !!---------------------------------------------------------------------- |
---|
117 | ! |
---|
118 | IF( ln_timing ) CALL timing_start('dia_wri') |
---|
119 | ! |
---|
120 | ! Output the initial state and forcings |
---|
121 | IF( ninist == 1 ) THEN |
---|
122 | CALL dia_wri_state( 'output.init' ) |
---|
123 | ninist = 0 |
---|
124 | ENDIF |
---|
125 | |
---|
126 | ! Output of initial vertical scale factor |
---|
127 | CALL iom_put("e3t_0", e3t_0(:,:,:) ) |
---|
128 | CALL iom_put("e3u_0", e3u_0(:,:,:) ) |
---|
129 | CALL iom_put("e3v_0", e3v_0(:,:,:) ) |
---|
130 | ! |
---|
131 | CALL iom_put( "e3t" , e3t_n(:,:,:) ) |
---|
132 | CALL iom_put( "e3u" , e3u_n(:,:,:) ) |
---|
133 | CALL iom_put( "e3v" , e3v_n(:,:,:) ) |
---|
134 | CALL iom_put( "e3w" , e3w_n(:,:,:) ) |
---|
135 | IF( iom_use("e3tdef") ) & |
---|
136 | CALL iom_put( "e3tdef" , ( ( e3t_n(:,:,:) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 ) |
---|
137 | |
---|
138 | IF( ll_wd ) THEN |
---|
139 | CALL iom_put( "ssh" , (sshn+ssh_ref)*tmask(:,:,1) ) ! sea surface height (brought back to the reference used for wetting and drying) |
---|
140 | ELSE |
---|
141 | CALL iom_put( "ssh" , sshn ) ! sea surface height |
---|
142 | ENDIF |
---|
143 | |
---|
144 | IF( iom_use("wetdep") ) & ! wet depth |
---|
145 | CALL iom_put( "wetdep" , ht_0(:,:) + sshn(:,:) ) |
---|
146 | |
---|
147 | CALL iom_put( "toce", tsn(:,:,:,jp_tem) ) ! 3D temperature |
---|
148 | CALL iom_put( "sst", tsn(:,:,1,jp_tem) ) ! surface temperature |
---|
149 | IF ( iom_use("sbt") ) THEN |
---|
150 | DO jj = 1, jpj |
---|
151 | DO ji = 1, jpi |
---|
152 | ikbot = mbkt(ji,jj) |
---|
153 | z2d(ji,jj) = tsn(ji,jj,ikbot,jp_tem) |
---|
154 | END DO |
---|
155 | END DO |
---|
156 | CALL iom_put( "sbt", z2d ) ! bottom temperature |
---|
157 | ENDIF |
---|
158 | |
---|
159 | CALL iom_put( "soce", tsn(:,:,:,jp_sal) ) ! 3D salinity |
---|
160 | CALL iom_put( "sss", tsn(:,:,1,jp_sal) ) ! surface salinity |
---|
161 | IF ( iom_use("sbs") ) THEN |
---|
162 | DO jj = 1, jpj |
---|
163 | DO ji = 1, jpi |
---|
164 | ikbot = mbkt(ji,jj) |
---|
165 | z2d(ji,jj) = tsn(ji,jj,ikbot,jp_sal) |
---|
166 | END DO |
---|
167 | END DO |
---|
168 | CALL iom_put( "sbs", z2d ) ! bottom salinity |
---|
169 | ENDIF |
---|
170 | |
---|
171 | IF ( iom_use("taubot") ) THEN ! bottom stress |
---|
172 | zztmp = rau0 * 0.25 |
---|
173 | z2d(:,:) = 0._wp |
---|
174 | DO jj = 2, jpjm1 |
---|
175 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
176 | zztmp2 = ( ( rCdU_bot(ji+1,jj)+rCdU_bot(ji ,jj) ) * un(ji ,jj,mbku(ji ,jj)) )**2 & |
---|
177 | & + ( ( rCdU_bot(ji ,jj)+rCdU_bot(ji-1,jj) ) * un(ji-1,jj,mbku(ji-1,jj)) )**2 & |
---|
178 | & + ( ( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj ) ) * vn(ji,jj ,mbkv(ji,jj )) )**2 & |
---|
179 | & + ( ( rCdU_bot(ji,jj )+rCdU_bot(ji,jj-1) ) * vn(ji,jj-1,mbkv(ji,jj-1)) )**2 |
---|
180 | z2d(ji,jj) = zztmp * SQRT( zztmp2 ) * tmask(ji,jj,1) |
---|
181 | ! |
---|
182 | END DO |
---|
183 | END DO |
---|
184 | CALL lbc_lnk( 'diawri', z2d, 'T', 1. ) |
---|
185 | CALL iom_put( "taubot", z2d ) |
---|
186 | ENDIF |
---|
187 | |
---|
188 | CALL iom_put( "uoce", un(:,:,:) ) ! 3D i-current |
---|
189 | CALL iom_put( "ssu", un(:,:,1) ) ! surface i-current |
---|
190 | IF ( iom_use("sbu") ) THEN |
---|
191 | DO jj = 1, jpj |
---|
192 | DO ji = 1, jpi |
---|
193 | ikbot = mbku(ji,jj) |
---|
194 | z2d(ji,jj) = un(ji,jj,ikbot) |
---|
195 | END DO |
---|
196 | END DO |
---|
197 | CALL iom_put( "sbu", z2d ) ! bottom i-current |
---|
198 | ENDIF |
---|
199 | |
---|
200 | CALL iom_put( "voce", vn(:,:,:) ) ! 3D j-current |
---|
201 | CALL iom_put( "ssv", vn(:,:,1) ) ! surface j-current |
---|
202 | IF ( iom_use("sbv") ) THEN |
---|
203 | DO jj = 1, jpj |
---|
204 | DO ji = 1, jpi |
---|
205 | ikbot = mbkv(ji,jj) |
---|
206 | z2d(ji,jj) = vn(ji,jj,ikbot) |
---|
207 | END DO |
---|
208 | END DO |
---|
209 | CALL iom_put( "sbv", z2d ) ! bottom j-current |
---|
210 | ENDIF |
---|
211 | |
---|
212 | CALL iom_put( "woce", wn ) ! vertical velocity |
---|
213 | IF( iom_use('w_masstr') .OR. iom_use('w_masstr2') ) THEN ! vertical mass transport & its square value |
---|
214 | ! Caution: in the VVL case, it only correponds to the baroclinic mass transport. |
---|
215 | z2d(:,:) = rau0 * e1e2t(:,:) |
---|
216 | DO jk = 1, jpk |
---|
217 | z3d(:,:,jk) = wn(:,:,jk) * z2d(:,:) |
---|
218 | END DO |
---|
219 | CALL iom_put( "w_masstr" , z3d ) |
---|
220 | IF( iom_use('w_masstr2') ) CALL iom_put( "w_masstr2", z3d(:,:,:) * z3d(:,:,:) ) |
---|
221 | ENDIF |
---|
222 | |
---|
223 | CALL iom_put( "avt" , avt ) ! T vert. eddy diff. coef. |
---|
224 | CALL iom_put( "avs" , avs ) ! S vert. eddy diff. coef. |
---|
225 | CALL iom_put( "avm" , avm ) ! T vert. eddy visc. coef. |
---|
226 | |
---|
227 | IF( iom_use('logavt') ) CALL iom_put( "logavt", LOG( MAX( 1.e-20_wp, avt(:,:,:) ) ) ) |
---|
228 | IF( iom_use('logavs') ) CALL iom_put( "logavs", LOG( MAX( 1.e-20_wp, avs(:,:,:) ) ) ) |
---|
229 | |
---|
230 | IF ( iom_use("sstgrad") .OR. iom_use("sstgrad2") ) THEN |
---|
231 | DO jj = 2, jpjm1 ! sst gradient |
---|
232 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
233 | zztmp = tsn(ji,jj,1,jp_tem) |
---|
234 | zztmpx = ( tsn(ji+1,jj,1,jp_tem) - zztmp ) * r1_e1u(ji,jj) + ( zztmp - tsn(ji-1,jj ,1,jp_tem) ) * r1_e1u(ji-1,jj) |
---|
235 | zztmpy = ( tsn(ji,jj+1,1,jp_tem) - zztmp ) * r1_e2v(ji,jj) + ( zztmp - tsn(ji ,jj-1,1,jp_tem) ) * r1_e2v(ji,jj-1) |
---|
236 | z2d(ji,jj) = 0.25 * ( zztmpx * zztmpx + zztmpy * zztmpy ) & |
---|
237 | & * umask(ji,jj,1) * umask(ji-1,jj,1) * vmask(ji,jj,1) * umask(ji,jj-1,1) |
---|
238 | END DO |
---|
239 | END DO |
---|
240 | CALL lbc_lnk( 'diawri', z2d, 'T', 1. ) |
---|
241 | CALL iom_put( "sstgrad2", z2d ) ! square of module of sst gradient |
---|
242 | z2d(:,:) = SQRT( z2d(:,:) ) |
---|
243 | CALL iom_put( "sstgrad" , z2d ) ! module of sst gradient |
---|
244 | ENDIF |
---|
245 | |
---|
246 | ! heat and salt contents |
---|
247 | IF( iom_use("heatc") ) THEN |
---|
248 | z2d(:,:) = 0._wp |
---|
249 | DO jk = 1, jpkm1 |
---|
250 | DO jj = 1, jpj |
---|
251 | DO ji = 1, jpi |
---|
252 | z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_tem) * tmask(ji,jj,jk) |
---|
253 | END DO |
---|
254 | END DO |
---|
255 | END DO |
---|
256 | CALL iom_put( "heatc", rau0_rcp * z2d ) ! vertically integrated heat content (J/m2) |
---|
257 | ENDIF |
---|
258 | |
---|
259 | IF( iom_use("saltc") ) THEN |
---|
260 | z2d(:,:) = 0._wp |
---|
261 | DO jk = 1, jpkm1 |
---|
262 | DO jj = 1, jpj |
---|
263 | DO ji = 1, jpi |
---|
264 | z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) * tmask(ji,jj,jk) |
---|
265 | END DO |
---|
266 | END DO |
---|
267 | END DO |
---|
268 | CALL iom_put( "saltc", rau0 * z2d ) ! vertically integrated salt content (PSU*kg/m2) |
---|
269 | ENDIF |
---|
270 | ! |
---|
271 | IF ( iom_use("eken") ) THEN |
---|
272 | z3d(:,:,jpk) = 0._wp |
---|
273 | DO jk = 1, jpkm1 |
---|
274 | DO jj = 2, jpjm1 |
---|
275 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
276 | zztmp = 0.25_wp * r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk) |
---|
277 | z3d(ji,jj,jk) = zztmp * ( un(ji-1,jj,jk)**2 * e2u(ji-1,jj) * e3u_n(ji-1,jj,jk) & |
---|
278 | & + un(ji ,jj,jk)**2 * e2u(ji ,jj) * e3u_n(ji ,jj,jk) & |
---|
279 | & + vn(ji,jj-1,jk)**2 * e1v(ji,jj-1) * e3v_n(ji,jj-1,jk) & |
---|
280 | & + vn(ji,jj ,jk)**2 * e1v(ji,jj ) * e3v_n(ji,jj ,jk) ) |
---|
281 | END DO |
---|
282 | END DO |
---|
283 | END DO |
---|
284 | CALL lbc_lnk( 'diawri', z3d, 'T', 1. ) |
---|
285 | CALL iom_put( "eken", z3d ) ! kinetic energy |
---|
286 | ENDIF |
---|
287 | ! |
---|
288 | CALL iom_put( "hdiv", hdivn ) ! Horizontal divergence |
---|
289 | ! |
---|
290 | IF( iom_use("u_masstr") .OR. iom_use("u_masstr_vint") .OR. iom_use("u_heattr") .OR. iom_use("u_salttr") ) THEN |
---|
291 | z3d(:,:,jpk) = 0.e0 |
---|
292 | z2d(:,:) = 0.e0 |
---|
293 | DO jk = 1, jpkm1 |
---|
294 | z3d(:,:,jk) = rau0 * un(:,:,jk) * e2u(:,:) * e3u_n(:,:,jk) * umask(:,:,jk) |
---|
295 | z2d(:,:) = z2d(:,:) + z3d(:,:,jk) |
---|
296 | END DO |
---|
297 | CALL iom_put( "u_masstr" , z3d ) ! mass transport in i-direction |
---|
298 | CALL iom_put( "u_masstr_vint", z2d ) ! mass transport in i-direction vertical sum |
---|
299 | ENDIF |
---|
300 | |
---|
301 | IF( iom_use("u_heattr") ) THEN |
---|
302 | z2d(:,:) = 0._wp |
---|
303 | DO jk = 1, jpkm1 |
---|
304 | DO jj = 2, jpjm1 |
---|
305 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
306 | z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_tem) + tsn(ji+1,jj,jk,jp_tem) ) |
---|
307 | END DO |
---|
308 | END DO |
---|
309 | END DO |
---|
310 | CALL lbc_lnk( 'diawri', z2d, 'U', -1. ) |
---|
311 | CALL iom_put( "u_heattr", 0.5*rcp * z2d ) ! heat transport in i-direction |
---|
312 | ENDIF |
---|
313 | |
---|
314 | IF( iom_use("u_salttr") ) THEN |
---|
315 | z2d(:,:) = 0.e0 |
---|
316 | DO jk = 1, jpkm1 |
---|
317 | DO jj = 2, jpjm1 |
---|
318 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
319 | z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_sal) + tsn(ji+1,jj,jk,jp_sal) ) |
---|
320 | END DO |
---|
321 | END DO |
---|
322 | END DO |
---|
323 | CALL lbc_lnk( 'diawri', z2d, 'U', -1. ) |
---|
324 | CALL iom_put( "u_salttr", 0.5 * z2d ) ! heat transport in i-direction |
---|
325 | ENDIF |
---|
326 | |
---|
327 | |
---|
328 | IF( iom_use("v_masstr") .OR. iom_use("v_heattr") .OR. iom_use("v_salttr") ) THEN |
---|
329 | z3d(:,:,jpk) = 0.e0 |
---|
330 | DO jk = 1, jpkm1 |
---|
331 | z3d(:,:,jk) = rau0 * vn(:,:,jk) * e1v(:,:) * e3v_n(:,:,jk) * vmask(:,:,jk) |
---|
332 | END DO |
---|
333 | CALL iom_put( "v_masstr", z3d ) ! mass transport in j-direction |
---|
334 | ENDIF |
---|
335 | |
---|
336 | IF( iom_use("v_heattr") ) THEN |
---|
337 | z2d(:,:) = 0.e0 |
---|
338 | DO jk = 1, jpkm1 |
---|
339 | DO jj = 2, jpjm1 |
---|
340 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
341 | z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_tem) + tsn(ji,jj+1,jk,jp_tem) ) |
---|
342 | END DO |
---|
343 | END DO |
---|
344 | END DO |
---|
345 | CALL lbc_lnk( 'diawri', z2d, 'V', -1. ) |
---|
346 | CALL iom_put( "v_heattr", 0.5*rcp * z2d ) ! heat transport in j-direction |
---|
347 | ENDIF |
---|
348 | |
---|
349 | IF( iom_use("v_salttr") ) THEN |
---|
350 | z2d(:,:) = 0._wp |
---|
351 | DO jk = 1, jpkm1 |
---|
352 | DO jj = 2, jpjm1 |
---|
353 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
354 | z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_sal) + tsn(ji,jj+1,jk,jp_sal) ) |
---|
355 | END DO |
---|
356 | END DO |
---|
357 | END DO |
---|
358 | CALL lbc_lnk( 'diawri', z2d, 'V', -1. ) |
---|
359 | CALL iom_put( "v_salttr", 0.5 * z2d ) ! heat transport in j-direction |
---|
360 | ENDIF |
---|
361 | |
---|
362 | IF( iom_use("tosmint") ) THEN |
---|
363 | z2d(:,:) = 0._wp |
---|
364 | DO jk = 1, jpkm1 |
---|
365 | DO jj = 2, jpjm1 |
---|
366 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
367 | z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_tem) |
---|
368 | END DO |
---|
369 | END DO |
---|
370 | END DO |
---|
371 | CALL lbc_lnk( 'diawri', z2d, 'T', -1. ) |
---|
372 | CALL iom_put( "tosmint", rau0 * z2d ) ! Vertical integral of temperature |
---|
373 | ENDIF |
---|
374 | IF( iom_use("somint") ) THEN |
---|
375 | z2d(:,:)=0._wp |
---|
376 | DO jk = 1, jpkm1 |
---|
377 | DO jj = 2, jpjm1 |
---|
378 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
379 | z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) |
---|
380 | END DO |
---|
381 | END DO |
---|
382 | END DO |
---|
383 | CALL lbc_lnk( 'diawri', z2d, 'T', -1. ) |
---|
384 | CALL iom_put( "somint", rau0 * z2d ) ! Vertical integral of salinity |
---|
385 | ENDIF |
---|
386 | |
---|
387 | CALL iom_put( "bn2", rn2 ) ! Brunt-Vaisala buoyancy frequency (N^2) |
---|
388 | ! |
---|
389 | |
---|
390 | IF (ln_diatmb) CALL dia_tmb ! tmb values |
---|
391 | |
---|
392 | IF (ln_dia25h) CALL dia_25h( kt ) ! 25h averaging |
---|
393 | |
---|
394 | IF( ln_timing ) CALL timing_stop('dia_wri') |
---|
395 | ! |
---|
396 | END SUBROUTINE dia_wri |
---|
397 | |
---|
398 | #else |
---|
399 | !!---------------------------------------------------------------------- |
---|
400 | !! Default option use IOIPSL library |
---|
401 | !!---------------------------------------------------------------------- |
---|
402 | |
---|
403 | INTEGER FUNCTION dia_wri_alloc() |
---|
404 | ! |
---|
405 | dia_wri_alloc = 0 |
---|
406 | ! |
---|
407 | END FUNCTION dia_wri_alloc |
---|
408 | |
---|
409 | SUBROUTINE dia_wri( kt ) |
---|
410 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
411 | |
---|
412 | IF( ninist == 1 ) THEN !== Output the initial state and forcings ==! |
---|
413 | CALL dia_wri_state( 'output.init' ) |
---|
414 | ninist = 0 |
---|
415 | ENDIF |
---|
416 | |
---|
417 | END SUBROUTINE dia_wri |
---|
418 | |
---|
419 | #endif |
---|
420 | |
---|
421 | SUBROUTINE dia_wri_state( cdfile_name ) |
---|
422 | !!--------------------------------------------------------------------- |
---|
423 | !! *** ROUTINE dia_wri_state *** |
---|
424 | !! |
---|
425 | !! ** Purpose : create a NetCDF file named cdfile_name which contains |
---|
426 | !! the instantaneous ocean state and forcing fields. |
---|
427 | !! Used to find errors in the initial state or save the last |
---|
428 | !! ocean state in case of abnormal end of a simulation |
---|
429 | !! |
---|
430 | !! ** Method : NetCDF files using ioipsl |
---|
431 | !! File 'output.init.nc' is created if ninist = 1 (namelist) |
---|
432 | !! File 'output.abort.nc' is created in case of abnormal job end |
---|
433 | !!---------------------------------------------------------------------- |
---|
434 | CHARACTER (len=* ), INTENT( in ) :: cdfile_name ! name of the file created |
---|
435 | !! |
---|
436 | INTEGER :: inum |
---|
437 | !!---------------------------------------------------------------------- |
---|
438 | ! |
---|
439 | IF(lwp) WRITE(numout,*) |
---|
440 | IF(lwp) WRITE(numout,*) 'dia_wri_state : single instantaneous ocean state' |
---|
441 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~ and forcing fields file created ' |
---|
442 | IF(lwp) WRITE(numout,*) ' and named :', cdfile_name, '...nc' |
---|
443 | |
---|
444 | #if defined key_si3 |
---|
445 | CALL iom_open( TRIM(cdfile_name), inum, ldwrt = .TRUE., kdlev = jpl ) |
---|
446 | #else |
---|
447 | CALL iom_open( TRIM(cdfile_name), inum, ldwrt = .TRUE. ) |
---|
448 | #endif |
---|
449 | |
---|
450 | CALL iom_rstput( 0, 0, inum, 'votemper', tsn(:,:,:,jp_tem) ) ! now temperature |
---|
451 | CALL iom_rstput( 0, 0, inum, 'vosaline', tsn(:,:,:,jp_sal) ) ! now salinity |
---|
452 | CALL iom_rstput( 0, 0, inum, 'sossheig', sshn ) ! sea surface height |
---|
453 | CALL iom_rstput( 0, 0, inum, 'vozocrtx', un ) ! now i-velocity |
---|
454 | CALL iom_rstput( 0, 0, inum, 'vomecrty', vn ) ! now j-velocity |
---|
455 | CALL iom_rstput( 0, 0, inum, 'vovecrtz', wn ) ! now k-velocity |
---|
456 | IF( ALLOCATED(ahtu) ) THEN |
---|
457 | CALL iom_rstput( 0, 0, inum, 'ahtu', ahtu ) ! aht at u-point |
---|
458 | CALL iom_rstput( 0, 0, inum, 'ahtv', ahtv ) ! aht at v-point |
---|
459 | ENDIF |
---|
460 | IF( ALLOCATED(ahmt) ) THEN |
---|
461 | CALL iom_rstput( 0, 0, inum, 'ahmt', ahmt ) ! ahmt at u-point |
---|
462 | CALL iom_rstput( 0, 0, inum, 'ahmf', ahmf ) ! ahmf at v-point |
---|
463 | ENDIF |
---|
464 | CALL iom_rstput( 0, 0, inum, 'sowaflup', emp - rnf ) ! freshwater budget |
---|
465 | CALL iom_rstput( 0, 0, inum, 'sohefldo', qsr + qns ) ! total heat flux |
---|
466 | CALL iom_rstput( 0, 0, inum, 'soshfldo', qsr ) ! solar heat flux |
---|
467 | CALL iom_rstput( 0, 0, inum, 'soicecov', fr_i ) ! ice fraction |
---|
468 | CALL iom_rstput( 0, 0, inum, 'sozotaux', utau ) ! i-wind stress |
---|
469 | CALL iom_rstput( 0, 0, inum, 'sometauy', vtau ) ! j-wind stress |
---|
470 | IF( .NOT.ln_linssh ) THEN |
---|
471 | CALL iom_rstput( 0, 0, inum, 'vovvldep', gdept_n ) ! T-cell depth |
---|
472 | CALL iom_rstput( 0, 0, inum, 'vovvle3t', e3t_n ) ! T-cell thickness |
---|
473 | END IF |
---|
474 | IF( ln_wave .AND. ln_sdw ) THEN |
---|
475 | CALL iom_rstput( 0, 0, inum, 'sdzocrtx', usd ) ! now StokesDrift i-velocity |
---|
476 | CALL iom_rstput( 0, 0, inum, 'sdmecrty', vsd ) ! now StokesDrift j-velocity |
---|
477 | CALL iom_rstput( 0, 0, inum, 'sdvecrtz', wsd ) ! now StokesDrift k-velocity |
---|
478 | ENDIF |
---|
479 | |
---|
480 | #if defined key_si3 |
---|
481 | IF( nn_ice == 2 ) THEN ! condition needed in case agrif + ice-model but no-ice in child grid |
---|
482 | CALL ice_wri_state( inum ) |
---|
483 | ENDIF |
---|
484 | #endif |
---|
485 | ! |
---|
486 | CALL iom_close( inum ) |
---|
487 | ! |
---|
488 | END SUBROUTINE dia_wri_state |
---|
489 | |
---|
490 | !!====================================================================== |
---|
491 | END MODULE diawri |
---|