1 | MODULE obs_write |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE obs_write *** |
---|
4 | !! Observation diagnosticss: Write observation related diagnostics |
---|
5 | !!===================================================================== |
---|
6 | |
---|
7 | !!---------------------------------------------------------------------- |
---|
8 | !! cdf_wri_p3d : Write profile observation diagnostics in NetCDF format |
---|
9 | !! obs_wri_sla : Write SLA observation related diagnostics |
---|
10 | !! obs_wri_sst : Write SST observation related diagnostics |
---|
11 | !! obs_wri_seaice: Write seaice observation related diagnostics |
---|
12 | !! cdf_wri_vel : Write velocity observation diagnostics in NetCDF format |
---|
13 | !!---------------------------------------------------------------------- |
---|
14 | |
---|
15 | !! * Modules used |
---|
16 | USE par_kind, ONLY : & ! Precision variables |
---|
17 | & wp |
---|
18 | USE in_out_manager ! I/O manager |
---|
19 | USE dom_oce ! Ocean space and time domain variables |
---|
20 | USE obs_types ! Observation type integer to character translation |
---|
21 | USE julian, ONLY : & ! Julian date routines |
---|
22 | & greg2jul |
---|
23 | USE obs_utils, ONLY : & ! Observation operator utility functions |
---|
24 | & chkerr |
---|
25 | USE obs_profiles_def ! Type definitions for profiles |
---|
26 | USE obs_surf_def ! Type defintions for surface observations |
---|
27 | USE obs_fbm ! Observation feedback I/O |
---|
28 | USE obs_grid ! Grid tools |
---|
29 | USE obs_conv ! Conversion between units |
---|
30 | USE obs_const |
---|
31 | USE obs_sla_types |
---|
32 | USE obs_rot_vel ! Rotation of velocities |
---|
33 | |
---|
34 | IMPLICIT NONE |
---|
35 | |
---|
36 | !! * Routine accessibility |
---|
37 | PRIVATE |
---|
38 | PUBLIC obs_wri_p3d, & ! Write profile observation related diagnostics |
---|
39 | & obs_wri_sla, & ! Write SLA observation related diagnostics |
---|
40 | & obs_wri_sst, & ! Write SST observation related diagnostics |
---|
41 | & obs_wri_sss, & ! Write SSS observation related diagnostics |
---|
42 | & obs_wri_seaice, & ! Write seaice observation related diagnostics |
---|
43 | & obs_wri_vel ! Write velocity observation related diagnostics |
---|
44 | |
---|
45 | CONTAINS |
---|
46 | |
---|
47 | SUBROUTINE obs_wri_p3d( cprefix, profdata ) |
---|
48 | !!----------------------------------------------------------------------- |
---|
49 | !! |
---|
50 | !! *** ROUTINE obs_wri_p3d *** |
---|
51 | !! |
---|
52 | !! ** Purpose : Write temperature and salinity (profile) observation |
---|
53 | !! related diagnostics |
---|
54 | !! |
---|
55 | !! ** Method : NetCDF |
---|
56 | !! |
---|
57 | !! ** Action : |
---|
58 | !! |
---|
59 | !! History : |
---|
60 | !! ! 06-04 (A. Vidard) Original |
---|
61 | !! ! 06-04 (A. Vidard) Reformatted |
---|
62 | !! ! 06-10 (A. Weaver) Cleanup |
---|
63 | !! ! 07-01 (K. Mogensen) Use profile data types |
---|
64 | !! ! 07-03 (K. Mogensen) General handling of profiles |
---|
65 | !! ! 09-01 (K. Mogensen) New feedback format |
---|
66 | !!----------------------------------------------------------------------- |
---|
67 | |
---|
68 | !! * Modules used |
---|
69 | |
---|
70 | !! * Arguments |
---|
71 | CHARACTER(LEN=*), INTENT(IN) :: & |
---|
72 | & cprefix ! Prefix for output files |
---|
73 | TYPE(obs_prof), INTENT(INOUT) :: & |
---|
74 | & profdata ! Full set of profile data |
---|
75 | |
---|
76 | !! * Local declarations |
---|
77 | TYPE(obfbdata) :: fbdata |
---|
78 | CHARACTER(LEN=40) :: & |
---|
79 | & cfname |
---|
80 | INTEGER :: & |
---|
81 | & ilevel |
---|
82 | INTEGER :: & |
---|
83 | & jvar, & |
---|
84 | & jo, & |
---|
85 | & jk, & |
---|
86 | & ik |
---|
87 | REAL(wp) :: & |
---|
88 | & zpres |
---|
89 | |
---|
90 | CALL init_obfbdata( fbdata ) |
---|
91 | |
---|
92 | ! Find maximum level |
---|
93 | ilevel = 0 |
---|
94 | DO jvar = 1, 2 |
---|
95 | ilevel = MAX( ilevel, MAXVAL( profdata%var(jvar)%nvlidx(:) ) ) |
---|
96 | ENDDO |
---|
97 | CALL alloc_obfbdata( fbdata, 2, profdata%nprof, ilevel, 1, 1, .TRUE. ) |
---|
98 | |
---|
99 | fbdata%cname(1) = 'POTM' |
---|
100 | fbdata%cname(2) = 'PSAL' |
---|
101 | fbdata%coblong(1) = 'Potential temperature' |
---|
102 | fbdata%coblong(2) = 'Practical salinity' |
---|
103 | fbdata%cobunit(1) = 'Degrees Celsius' |
---|
104 | fbdata%cobunit(2) = 'PSU' |
---|
105 | fbdata%cextname(1) = 'TEMP' |
---|
106 | fbdata%cextlong(1) = 'Insitu temperature' |
---|
107 | fbdata%cextunit(1) = 'Degrees Celsius' |
---|
108 | fbdata%caddname(1) = 'Hx' |
---|
109 | fbdata%caddlong(1,1) = 'Model interpolated potential temperature' |
---|
110 | fbdata%caddlong(1,2) = 'Model interpolated practical salinity' |
---|
111 | fbdata%caddunit(1,1) = 'Degrees Celsius' |
---|
112 | fbdata%caddunit(1,2) = 'PSU' |
---|
113 | |
---|
114 | WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc |
---|
115 | |
---|
116 | IF(lwp) THEN |
---|
117 | WRITE(numout,*) |
---|
118 | WRITE(numout,*)'obs_wri_p3d :' |
---|
119 | WRITE(numout,*)'~~~~~~~~~~~~~' |
---|
120 | WRITE(numout,*)'Writing profile feedback file : ',TRIM(cfname) |
---|
121 | ENDIF |
---|
122 | |
---|
123 | ! Transform obs_prof data structure into obfbdata structure |
---|
124 | fbdata%cdjuldref = '19500101000000' |
---|
125 | DO jo = 1, profdata%nprof |
---|
126 | fbdata%plam(jo) = profdata%rlam(jo) |
---|
127 | fbdata%pphi(jo) = profdata%rphi(jo) |
---|
128 | WRITE(fbdata%cdtyp(jo),'(I4)') profdata%ntyp(jo) |
---|
129 | fbdata%ivqc(jo,:) = profdata%ivqc(jo,:) |
---|
130 | fbdata%ivqcf(:,jo,:) = profdata%ivqcf(:,jo,:) |
---|
131 | IF ( profdata%nqc(jo) > 10 ) THEN |
---|
132 | fbdata%ioqc(jo) = 4 |
---|
133 | fbdata%ioqcf(1,jo) = profdata%nqcf(1,jo) |
---|
134 | fbdata%ioqcf(2,jo) = profdata%nqc(jo) - 10 |
---|
135 | ELSE |
---|
136 | fbdata%ioqc(jo) = profdata%nqc(jo) |
---|
137 | fbdata%ioqcf(:,jo) = profdata%nqcf(:,jo) |
---|
138 | ENDIF |
---|
139 | fbdata%ipqc(jo) = profdata%ipqc(jo) |
---|
140 | fbdata%ipqcf(:,jo) = profdata%ipqcf(:,jo) |
---|
141 | fbdata%itqc(jo) = profdata%itqc(jo) |
---|
142 | fbdata%itqcf(:,jo) = profdata%itqcf(:,jo) |
---|
143 | fbdata%cdwmo(jo) = profdata%cwmo(jo) |
---|
144 | fbdata%kindex(jo) = profdata%npfil(jo) |
---|
145 | DO jvar = 1, profdata%nvar |
---|
146 | IF (ln_grid_global) THEN |
---|
147 | fbdata%iobsi(jo,jvar) = profdata%mi(jo,jvar) |
---|
148 | fbdata%iobsj(jo,jvar) = profdata%mj(jo,jvar) |
---|
149 | ELSE |
---|
150 | fbdata%iobsi(jo,jvar) = mig(profdata%mi(jo,jvar)) |
---|
151 | fbdata%iobsj(jo,jvar) = mjg(profdata%mj(jo,jvar)) |
---|
152 | ENDIF |
---|
153 | ENDDO |
---|
154 | CALL greg2jul( 0, & |
---|
155 | & profdata%nmin(jo), & |
---|
156 | & profdata%nhou(jo), & |
---|
157 | & profdata%nday(jo), & |
---|
158 | & profdata%nmon(jo), & |
---|
159 | & profdata%nyea(jo), & |
---|
160 | & fbdata%ptim(jo), & |
---|
161 | & krefdate = 19500101 ) |
---|
162 | ! Reform the profiles arrays for output |
---|
163 | DO jvar = 1, 2 |
---|
164 | DO jk = profdata%npvsta(jo,jvar), profdata%npvend(jo,jvar) |
---|
165 | ik = profdata%var(jvar)%nvlidx(jk) |
---|
166 | fbdata%padd(ik,jo,1,jvar) = profdata%var(jvar)%vmod(jk) |
---|
167 | fbdata%pob(ik,jo,jvar) = profdata%var(jvar)%vobs(jk) |
---|
168 | fbdata%pdep(ik,jo) = profdata%var(jvar)%vdep(jk) |
---|
169 | fbdata%idqc(ik,jo) = profdata%var(jvar)%idqc(jk) |
---|
170 | fbdata%idqcf(:,ik,jo) = profdata%var(jvar)%idqcf(:,jk) |
---|
171 | IF ( profdata%var(jvar)%nvqc(jk) > 10 ) THEN |
---|
172 | fbdata%ivlqc(ik,jo,jvar) = 4 |
---|
173 | fbdata%ivlqcf(1,ik,jo,jvar) = profdata%var(jvar)%nvqcf(1,jk) |
---|
174 | fbdata%ivlqcf(2,ik,jo,jvar) = profdata%var(jvar)%nvqc(jk) - 10 |
---|
175 | ELSE |
---|
176 | fbdata%ivlqc(ik,jo,jvar) = profdata%var(jvar)%nvqc(jk) |
---|
177 | fbdata%ivlqcf(:,ik,jo,jvar) = profdata%var(jvar)%nvqcf(:,jk) |
---|
178 | ENDIF |
---|
179 | fbdata%iobsk(ik,jo,jvar) = profdata%var(jvar)%mvk(jk) |
---|
180 | IF ( jvar == 1 ) THEN |
---|
181 | fbdata%pext(ik,jo,1) = profdata%var(jvar)%vext(jk,1) |
---|
182 | ENDIF |
---|
183 | END DO |
---|
184 | END DO |
---|
185 | END DO |
---|
186 | |
---|
187 | ! Convert insitu temperature to potential temperature using the model |
---|
188 | ! salinity if no potential temperature |
---|
189 | DO jo = 1, fbdata%nobs |
---|
190 | IF ( fbdata%pphi(jo) < 9999.0 ) THEN |
---|
191 | DO jk = 1, fbdata%nlev |
---|
192 | IF ( ( fbdata%pob(jk,jo,1) >= 9999.0 ) .AND. & |
---|
193 | & ( fbdata%pdep(jk,jo) < 9999.0 ) .AND. & |
---|
194 | & ( fbdata%padd(jk,jo,1,2) < 9999.0 ) .AND. & |
---|
195 | & ( fbdata%pext(jk,jo,1) < 9999.0 ) ) THEN |
---|
196 | zpres = dep_to_p( REAL(fbdata%pdep(jk,jo),wp), & |
---|
197 | & REAL(fbdata%pphi(jo),wp) ) |
---|
198 | fbdata%pob(jk,jo,1) = potemp( & |
---|
199 | & REAL(fbdata%padd(jk,jo,1,2), wp), & |
---|
200 | & REAL(fbdata%pext(jk,jo,1), wp), & |
---|
201 | & zpres, 0.0_wp ) |
---|
202 | ENDIF |
---|
203 | ENDDO |
---|
204 | ENDIF |
---|
205 | ENDDO |
---|
206 | |
---|
207 | ! Write the obfbdata structure |
---|
208 | CALL write_obfbdata( cfname, fbdata ) |
---|
209 | |
---|
210 | CALL dealloc_obfbdata( fbdata ) |
---|
211 | |
---|
212 | END SUBROUTINE obs_wri_p3d |
---|
213 | |
---|
214 | SUBROUTINE obs_wri_sla( cprefix, sladata ) |
---|
215 | !!----------------------------------------------------------------------- |
---|
216 | !! |
---|
217 | !! *** ROUTINE obs_wri_sla *** |
---|
218 | !! |
---|
219 | !! ** Purpose : Write SLA observation diagnostics |
---|
220 | !! related |
---|
221 | !! |
---|
222 | !! ** Method : NetCDF |
---|
223 | !! |
---|
224 | !! ** Action : |
---|
225 | !! |
---|
226 | !! ! 07-03 (K. Mogensen) Original |
---|
227 | !! ! 09-01 (K. Mogensen) New feedback format. |
---|
228 | !!----------------------------------------------------------------------- |
---|
229 | |
---|
230 | !! * Modules used |
---|
231 | IMPLICIT NONE |
---|
232 | |
---|
233 | !! * Arguments |
---|
234 | CHARACTER(LEN=*), INTENT(IN) :: & |
---|
235 | & cprefix ! Prefix for output files |
---|
236 | TYPE(obs_surf), INTENT(INOUT) :: & |
---|
237 | & sladata ! Full set of SLAa |
---|
238 | |
---|
239 | !! * Local declarations |
---|
240 | TYPE(obfbdata) :: fbdata |
---|
241 | CHARACTER(LEN=40) :: & |
---|
242 | & cfname ! netCDF filename |
---|
243 | CHARACTER(LEN=12), PARAMETER :: & |
---|
244 | & cpname = 'obs_wri_sla' |
---|
245 | INTEGER :: & |
---|
246 | & jo |
---|
247 | |
---|
248 | CALL init_obfbdata( fbdata ) |
---|
249 | |
---|
250 | CALL alloc_obfbdata( fbdata, 1, sladata%nsurf, 1, 1, 2, .TRUE. ) |
---|
251 | |
---|
252 | fbdata%cname(1) = 'SLA' |
---|
253 | fbdata%coblong(1) = 'Sea level anomaly' |
---|
254 | fbdata%cobunit(1) = 'metre' |
---|
255 | fbdata%cextname(1) = 'SSH' |
---|
256 | fbdata%cextlong(1) = 'Model Sea surface height' |
---|
257 | fbdata%cextunit(1) = 'metre' |
---|
258 | fbdata%cextname(2) = 'MDT' |
---|
259 | fbdata%cextlong(2) = 'Mean dynamic topography' |
---|
260 | fbdata%cextunit(2) = 'metre' |
---|
261 | fbdata%caddname(1) = 'Hx' |
---|
262 | fbdata%caddlong(1,1) = 'Model interpolated SSH - MDT' |
---|
263 | fbdata%caddunit(1,1) = 'metre' |
---|
264 | fbdata%cgrid(1) = 'T' |
---|
265 | |
---|
266 | WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc |
---|
267 | |
---|
268 | IF(lwp) THEN |
---|
269 | WRITE(numout,*) |
---|
270 | WRITE(numout,*)'obs_wri_sla :' |
---|
271 | WRITE(numout,*)'~~~~~~~~~~~~~' |
---|
272 | WRITE(numout,*)'Writing SLA feedback file : ',TRIM(cfname) |
---|
273 | ENDIF |
---|
274 | |
---|
275 | ! Transform obs_prof data structure into obfbdata structure |
---|
276 | fbdata%cdjuldref = '19500101000000' |
---|
277 | DO jo = 1, sladata%nsurf |
---|
278 | fbdata%plam(jo) = sladata%rlam(jo) |
---|
279 | fbdata%pphi(jo) = sladata%rphi(jo) |
---|
280 | WRITE(fbdata%cdtyp(jo),'(I4)') sladata%ntyp(jo) |
---|
281 | fbdata%ivqc(jo,:) = 0 |
---|
282 | fbdata%ivqcf(:,jo,:) = 0 |
---|
283 | IF ( sladata%nqc(jo) > 10 ) THEN |
---|
284 | fbdata%ioqc(jo) = 4 |
---|
285 | fbdata%ioqcf(1,jo) = 0 |
---|
286 | fbdata%ioqcf(2,jo) = sladata%nqc(jo) - 10 |
---|
287 | ELSE |
---|
288 | fbdata%ioqc(jo) = sladata%nqc(jo) |
---|
289 | fbdata%ioqcf(:,jo) = 0 |
---|
290 | ENDIF |
---|
291 | fbdata%ipqc(jo) = 0 |
---|
292 | fbdata%ipqcf(:,jo) = 0 |
---|
293 | fbdata%itqc(jo) = 0 |
---|
294 | fbdata%itqcf(:,jo) = 0 |
---|
295 | fbdata%cdwmo(jo) = cmissions(sladata%ntyp(jo)) |
---|
296 | fbdata%kindex(jo) = sladata%nsfil(jo) |
---|
297 | IF (ln_grid_global) THEN |
---|
298 | fbdata%iobsi(jo,1) = sladata%mi(jo) |
---|
299 | fbdata%iobsj(jo,1) = sladata%mj(jo) |
---|
300 | ELSE |
---|
301 | fbdata%iobsi(jo,1) = mig(sladata%mi(jo)) |
---|
302 | fbdata%iobsj(jo,1) = mjg(sladata%mj(jo)) |
---|
303 | ENDIF |
---|
304 | CALL greg2jul( 0, & |
---|
305 | & sladata%nmin(jo), & |
---|
306 | & sladata%nhou(jo), & |
---|
307 | & sladata%nday(jo), & |
---|
308 | & sladata%nmon(jo), & |
---|
309 | & sladata%nyea(jo), & |
---|
310 | & fbdata%ptim(jo), & |
---|
311 | & krefdate = 19500101 ) |
---|
312 | fbdata%padd(1,jo,1,1) = sladata%rmod(jo,1) |
---|
313 | fbdata%pob(1,jo,1) = sladata%robs(jo,1) |
---|
314 | fbdata%pdep(1,jo) = 0.0 |
---|
315 | fbdata%idqc(1,jo) = 0 |
---|
316 | fbdata%idqcf(:,1,jo) = 0 |
---|
317 | IF ( sladata%nqc(jo) > 10 ) THEN |
---|
318 | fbdata%ivlqc(1,jo,1) = 4 |
---|
319 | fbdata%ivlqcf(1,1,jo,1) = 0 |
---|
320 | fbdata%ivlqcf(2,1,jo,1) = sladata%nqc(jo) - 10 |
---|
321 | ELSE |
---|
322 | fbdata%ivlqc(1,jo,1) = sladata%nqc(jo) |
---|
323 | fbdata%ivlqcf(:,1,jo,1) = 0 |
---|
324 | ENDIF |
---|
325 | fbdata%iobsk(1,jo,1) = 0 |
---|
326 | fbdata%pext(1,jo,1) = sladata%rext(jo,1) |
---|
327 | fbdata%pext(1,jo,2) = sladata%rext(jo,2) |
---|
328 | |
---|
329 | END DO |
---|
330 | |
---|
331 | ! Write the obfbdata structure |
---|
332 | CALL write_obfbdata( cfname, fbdata ) |
---|
333 | |
---|
334 | CALL dealloc_obfbdata( fbdata ) |
---|
335 | |
---|
336 | END SUBROUTINE obs_wri_sla |
---|
337 | |
---|
338 | SUBROUTINE obs_wri_sst( cprefix, sstdata ) |
---|
339 | !!----------------------------------------------------------------------- |
---|
340 | !! |
---|
341 | !! *** ROUTINE obs_wri_sst *** |
---|
342 | !! |
---|
343 | !! ** Purpose : Write SST observation diagnostics |
---|
344 | !! related |
---|
345 | !! |
---|
346 | !! ** Method : NetCDF |
---|
347 | !! |
---|
348 | !! ** Action : |
---|
349 | !! |
---|
350 | !! ! 07-07 (S. Ricci) Original |
---|
351 | !! ! 09-01 (K. Mogensen) New feedback format. |
---|
352 | !!----------------------------------------------------------------------- |
---|
353 | |
---|
354 | !! * Modules used |
---|
355 | IMPLICIT NONE |
---|
356 | |
---|
357 | !! * Arguments |
---|
358 | CHARACTER(LEN=*), INTENT(IN) :: & |
---|
359 | & cprefix ! Prefix for output files |
---|
360 | TYPE(obs_surf), INTENT(INOUT) :: & |
---|
361 | & sstdata ! Full set of SST |
---|
362 | |
---|
363 | !! * Local declarations |
---|
364 | TYPE(obfbdata) :: fbdata |
---|
365 | CHARACTER(LEN=40) :: & |
---|
366 | & cfname ! netCDF filename |
---|
367 | CHARACTER(LEN=12), PARAMETER :: & |
---|
368 | & cpname = 'obs_wri_sst' |
---|
369 | INTEGER :: & |
---|
370 | & jo |
---|
371 | |
---|
372 | CALL init_obfbdata( fbdata ) |
---|
373 | |
---|
374 | CALL alloc_obfbdata( fbdata, 1, sstdata%nsurf, 1, 1, 0, .TRUE. ) |
---|
375 | |
---|
376 | fbdata%cname(1) = 'SST' |
---|
377 | fbdata%coblong(1) = 'Sea surface temperature' |
---|
378 | fbdata%cobunit(1) = 'Degree centigrade' |
---|
379 | fbdata%caddname(1) = 'Hx' |
---|
380 | fbdata%caddlong(1,1) = 'Model interpolated SST' |
---|
381 | fbdata%caddunit(1,1) = 'Degree centigrade' |
---|
382 | |
---|
383 | WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc |
---|
384 | |
---|
385 | IF(lwp) THEN |
---|
386 | WRITE(numout,*) |
---|
387 | WRITE(numout,*)'obs_wri_sst :' |
---|
388 | WRITE(numout,*)'~~~~~~~~~~~~~' |
---|
389 | WRITE(numout,*)'Writing SST feedback file : ',TRIM(cfname) |
---|
390 | ENDIF |
---|
391 | |
---|
392 | ! Transform obs_prof data structure into obfbdata structure |
---|
393 | fbdata%cdjuldref = '19500101000000' |
---|
394 | DO jo = 1, sstdata%nsurf |
---|
395 | fbdata%plam(jo) = sstdata%rlam(jo) |
---|
396 | fbdata%pphi(jo) = sstdata%rphi(jo) |
---|
397 | WRITE(fbdata%cdtyp(jo),'(I4)') sstdata%ntyp(jo) |
---|
398 | fbdata%ivqc(jo,:) = 0 |
---|
399 | fbdata%ivqcf(:,jo,:) = 0 |
---|
400 | IF ( sstdata%nqc(jo) > 10 ) THEN |
---|
401 | fbdata%ioqc(jo) = 4 |
---|
402 | fbdata%ioqcf(1,jo) = 0 |
---|
403 | fbdata%ioqcf(2,jo) = sstdata%nqc(jo) - 10 |
---|
404 | ELSE |
---|
405 | fbdata%ioqc(jo) = MAX(sstdata%nqc(jo),1) |
---|
406 | fbdata%ioqcf(:,jo) = 0 |
---|
407 | ENDIF |
---|
408 | fbdata%ipqc(jo) = 0 |
---|
409 | fbdata%ipqcf(:,jo) = 0 |
---|
410 | fbdata%itqc(jo) = 0 |
---|
411 | fbdata%itqcf(:,jo) = 0 |
---|
412 | fbdata%cdwmo(jo) = '' |
---|
413 | fbdata%kindex(jo) = sstdata%nsfil(jo) |
---|
414 | IF (ln_grid_global) THEN |
---|
415 | fbdata%iobsi(jo,1) = sstdata%mi(jo) |
---|
416 | fbdata%iobsj(jo,1) = sstdata%mj(jo) |
---|
417 | ELSE |
---|
418 | fbdata%iobsi(jo,1) = mig(sstdata%mi(jo)) |
---|
419 | fbdata%iobsj(jo,1) = mjg(sstdata%mj(jo)) |
---|
420 | ENDIF |
---|
421 | CALL greg2jul( 0, & |
---|
422 | & sstdata%nmin(jo), & |
---|
423 | & sstdata%nhou(jo), & |
---|
424 | & sstdata%nday(jo), & |
---|
425 | & sstdata%nmon(jo), & |
---|
426 | & sstdata%nyea(jo), & |
---|
427 | & fbdata%ptim(jo), & |
---|
428 | & krefdate = 19500101 ) |
---|
429 | fbdata%padd(1,jo,1,1) = sstdata%rmod(jo,1) |
---|
430 | fbdata%pob(1,jo,1) = sstdata%robs(jo,1) |
---|
431 | fbdata%pdep(1,jo) = 0.0 |
---|
432 | fbdata%idqc(1,jo) = 0 |
---|
433 | fbdata%idqcf(:,1,jo) = 0 |
---|
434 | IF ( sstdata%nqc(jo) > 10 ) THEN |
---|
435 | fbdata%ivlqc(1,jo,1) = 4 |
---|
436 | fbdata%ivlqcf(1,1,jo,1) = 0 |
---|
437 | fbdata%ivlqcf(2,1,jo,1) = sstdata%nqc(jo) - 10 |
---|
438 | ELSE |
---|
439 | fbdata%ivlqc(1,jo,1) = MAX(sstdata%nqc(jo),1) |
---|
440 | fbdata%ivlqcf(:,1,jo,1) = 0 |
---|
441 | ENDIF |
---|
442 | fbdata%iobsk(1,jo,1) = 0 |
---|
443 | |
---|
444 | END DO |
---|
445 | |
---|
446 | ! Write the obfbdata structure |
---|
447 | CALL write_obfbdata( cfname, fbdata ) |
---|
448 | |
---|
449 | CALL dealloc_obfbdata( fbdata ) |
---|
450 | |
---|
451 | END SUBROUTINE obs_wri_sst |
---|
452 | |
---|
453 | SUBROUTINE obs_wri_sss |
---|
454 | END SUBROUTINE obs_wri_sss |
---|
455 | |
---|
456 | SUBROUTINE obs_wri_seaice( cprefix, seaicedata ) |
---|
457 | !!----------------------------------------------------------------------- |
---|
458 | !! |
---|
459 | !! *** ROUTINE obs_wri_seaice *** |
---|
460 | !! |
---|
461 | !! ** Purpose : Write sea ice observation diagnostics |
---|
462 | !! related |
---|
463 | !! |
---|
464 | !! ** Method : NetCDF |
---|
465 | !! |
---|
466 | !! ** Action : |
---|
467 | !! |
---|
468 | !! ! 07-07 (S. Ricci) Original |
---|
469 | !! ! 09-01 (K. Mogensen) New feedback format. |
---|
470 | !!----------------------------------------------------------------------- |
---|
471 | |
---|
472 | !! * Modules used |
---|
473 | IMPLICIT NONE |
---|
474 | |
---|
475 | !! * Arguments |
---|
476 | CHARACTER(LEN=*), INTENT(IN) :: & |
---|
477 | & cprefix ! Prefix for output files |
---|
478 | TYPE(obs_surf), INTENT(INOUT) :: & |
---|
479 | & seaicedata ! Full set of sea ice |
---|
480 | TYPE(obfbdata) :: fbdata |
---|
481 | CHARACTER(LEN=40) :: & |
---|
482 | & cfname ! netCDF filename |
---|
483 | CHARACTER(LEN=12), PARAMETER :: & |
---|
484 | & cpname = 'obs_wri_seaice' |
---|
485 | INTEGER :: & |
---|
486 | & jo |
---|
487 | |
---|
488 | CALL init_obfbdata( fbdata ) |
---|
489 | |
---|
490 | CALL alloc_obfbdata( fbdata, 1, seaicedata%nsurf, 1, 1, 0, .TRUE. ) |
---|
491 | |
---|
492 | fbdata%cname(1) = 'SEAICE' |
---|
493 | fbdata%coblong(1) = 'Sea ice' |
---|
494 | fbdata%cobunit(1) = 'Fraction' |
---|
495 | fbdata%caddname(1) = 'Hx' |
---|
496 | fbdata%caddlong(1,1) = 'Model interpolated ICE' |
---|
497 | fbdata%caddunit(1,1) = 'Fraction' |
---|
498 | |
---|
499 | WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc |
---|
500 | |
---|
501 | IF(lwp) THEN |
---|
502 | WRITE(numout,*) |
---|
503 | WRITE(numout,*)'obs_wri_seaice :' |
---|
504 | WRITE(numout,*)'~~~~~~~~~~~~~~~~' |
---|
505 | WRITE(numout,*)'Writing SEAICE feedback file : ',TRIM(cfname) |
---|
506 | ENDIF |
---|
507 | |
---|
508 | ! Transform obs_prof data structure into obfbdata structure |
---|
509 | fbdata%cdjuldref = '19500101000000' |
---|
510 | DO jo = 1, seaicedata%nsurf |
---|
511 | fbdata%plam(jo) = seaicedata%rlam(jo) |
---|
512 | fbdata%pphi(jo) = seaicedata%rphi(jo) |
---|
513 | WRITE(fbdata%cdtyp(jo),'(I4)') seaicedata%ntyp(jo) |
---|
514 | fbdata%ivqc(jo,:) = 0 |
---|
515 | fbdata%ivqcf(:,jo,:) = 0 |
---|
516 | IF ( seaicedata%nqc(jo) > 10 ) THEN |
---|
517 | fbdata%ioqc(jo) = 4 |
---|
518 | fbdata%ioqcf(1,jo) = 0 |
---|
519 | fbdata%ioqcf(2,jo) = seaicedata%nqc(jo) - 10 |
---|
520 | ELSE |
---|
521 | fbdata%ioqc(jo) = MAX(seaicedata%nqc(jo),1) |
---|
522 | fbdata%ioqcf(:,jo) = 0 |
---|
523 | ENDIF |
---|
524 | fbdata%ipqc(jo) = 0 |
---|
525 | fbdata%ipqcf(:,jo) = 0 |
---|
526 | fbdata%itqc(jo) = 0 |
---|
527 | fbdata%itqcf(:,jo) = 0 |
---|
528 | fbdata%cdwmo(jo) = '' |
---|
529 | fbdata%kindex(jo) = seaicedata%nsfil(jo) |
---|
530 | IF (ln_grid_global) THEN |
---|
531 | fbdata%iobsi(jo,1) = seaicedata%mi(jo) |
---|
532 | fbdata%iobsj(jo,1) = seaicedata%mj(jo) |
---|
533 | ELSE |
---|
534 | fbdata%iobsi(jo,1) = mig(seaicedata%mi(jo)) |
---|
535 | fbdata%iobsj(jo,1) = mjg(seaicedata%mj(jo)) |
---|
536 | ENDIF |
---|
537 | CALL greg2jul( 0, & |
---|
538 | & seaicedata%nmin(jo), & |
---|
539 | & seaicedata%nhou(jo), & |
---|
540 | & seaicedata%nday(jo), & |
---|
541 | & seaicedata%nmon(jo), & |
---|
542 | & seaicedata%nyea(jo), & |
---|
543 | & fbdata%ptim(jo), & |
---|
544 | & krefdate = 19500101 ) |
---|
545 | fbdata%padd(1,jo,1,1) = seaicedata%rmod(jo,1) |
---|
546 | fbdata%pob(1,jo,1) = seaicedata%robs(jo,1) |
---|
547 | fbdata%pdep(1,jo) = 0.0 |
---|
548 | fbdata%idqc(1,jo) = 0 |
---|
549 | fbdata%idqcf(:,1,jo) = 0 |
---|
550 | IF ( seaicedata%nqc(jo) > 10 ) THEN |
---|
551 | fbdata%ivlqc(1,jo,1) = 4 |
---|
552 | fbdata%ivlqcf(1,1,jo,1) = 0 |
---|
553 | fbdata%ivlqcf(2,1,jo,1) = seaicedata%nqc(jo) - 10 |
---|
554 | ELSE |
---|
555 | fbdata%ivlqc(1,jo,1) = MAX(seaicedata%nqc(jo),1) |
---|
556 | fbdata%ivlqcf(:,1,jo,1) = 0 |
---|
557 | ENDIF |
---|
558 | fbdata%iobsk(1,jo,1) = 0 |
---|
559 | |
---|
560 | END DO |
---|
561 | |
---|
562 | ! Write the obfbdata structure |
---|
563 | CALL write_obfbdata( cfname, fbdata ) |
---|
564 | |
---|
565 | CALL dealloc_obfbdata( fbdata ) |
---|
566 | |
---|
567 | |
---|
568 | !! * Local declarations |
---|
569 | END SUBROUTINE obs_wri_seaice |
---|
570 | |
---|
571 | SUBROUTINE obs_wri_vel( cprefix, profdata, k2dint ) |
---|
572 | !!----------------------------------------------------------------------- |
---|
573 | !! |
---|
574 | !! *** ROUTINE obs_wri_p3d *** |
---|
575 | !! |
---|
576 | !! ** Purpose : Write current (profile) observation |
---|
577 | !! related diagnostics |
---|
578 | !! |
---|
579 | !! ** Method : NetCDF |
---|
580 | !! |
---|
581 | !! ** Action : |
---|
582 | !! |
---|
583 | !! History : |
---|
584 | !! ! 09-01 (K. Mogensen) New feedback format routine |
---|
585 | !!----------------------------------------------------------------------- |
---|
586 | |
---|
587 | !! * Modules used |
---|
588 | |
---|
589 | !! * Arguments |
---|
590 | CHARACTER(LEN=*), INTENT(IN) :: & |
---|
591 | & cprefix ! Prefix for output files |
---|
592 | TYPE(obs_prof), INTENT(INOUT) :: & |
---|
593 | & profdata ! Full set of profile data |
---|
594 | INTEGER, INTENT(IN) :: & |
---|
595 | & k2dint ! Horizontal interpolation method |
---|
596 | |
---|
597 | !! * Local declarations |
---|
598 | TYPE(obfbdata) :: fbdata |
---|
599 | CHARACTER(LEN=40) :: & |
---|
600 | & cfname |
---|
601 | INTEGER :: & |
---|
602 | & ilevel |
---|
603 | INTEGER :: & |
---|
604 | & jvar, & |
---|
605 | & jo, & |
---|
606 | & jk, & |
---|
607 | & ik |
---|
608 | REAL(wp) :: & |
---|
609 | & zpres |
---|
610 | REAL(wp), DIMENSION(:), ALLOCATABLE :: & |
---|
611 | & zu, & |
---|
612 | & zv |
---|
613 | |
---|
614 | CALL init_obfbdata( fbdata ) |
---|
615 | |
---|
616 | ! Find maximum level |
---|
617 | ilevel = 0 |
---|
618 | DO jvar = 1, 2 |
---|
619 | ilevel = MAX( ilevel, MAXVAL( profdata%var(jvar)%nvlidx(:) ) ) |
---|
620 | ENDDO |
---|
621 | CALL alloc_obfbdata( fbdata, 2, profdata%nprof, ilevel, 2, 0, .TRUE. ) |
---|
622 | |
---|
623 | fbdata%cname(1) = 'UVEL' |
---|
624 | fbdata%cname(2) = 'VVEL' |
---|
625 | fbdata%coblong(1) = 'Zonal velocity' |
---|
626 | fbdata%coblong(2) = 'Meridional velocity' |
---|
627 | fbdata%cobunit(1) = 'm/s' |
---|
628 | fbdata%cobunit(2) = 'm/s' |
---|
629 | fbdata%caddname(1) = 'Hx' |
---|
630 | fbdata%caddlong(1,1) = 'Model interpolated zonal velocity' |
---|
631 | fbdata%caddlong(1,2) = 'Model interpolated meridional velocity' |
---|
632 | fbdata%caddunit(1,1) = 'm/s' |
---|
633 | fbdata%caddunit(1,2) = 'm/s' |
---|
634 | fbdata%caddname(2) = 'HxGRID' |
---|
635 | fbdata%caddlong(2,1) = 'Model interpolated zonal velocity (model grid)' |
---|
636 | fbdata%caddlong(2,2) = 'Model interpolated meridional velocity (model grid)' |
---|
637 | fbdata%caddunit(2,1) = 'm/s' |
---|
638 | fbdata%caddunit(2,2) = 'm/s' |
---|
639 | |
---|
640 | WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc |
---|
641 | |
---|
642 | IF(lwp) THEN |
---|
643 | WRITE(numout,*) |
---|
644 | WRITE(numout,*)'obs_wri_vel :' |
---|
645 | WRITE(numout,*)'~~~~~~~~~~~~~' |
---|
646 | WRITE(numout,*)'Writing velocuty feedback file : ',TRIM(cfname) |
---|
647 | ENDIF |
---|
648 | |
---|
649 | ALLOCATE( & |
---|
650 | & zu(profdata%nvprot(1)), & |
---|
651 | & zv(profdata%nvprot(2)) & |
---|
652 | & ) |
---|
653 | CALL obs_rotvel( profdata, k2dint, zu, zv ) |
---|
654 | |
---|
655 | ! Transform obs_prof data structure into obfbdata structure |
---|
656 | fbdata%cdjuldref = '19500101000000' |
---|
657 | DO jo = 1, profdata%nprof |
---|
658 | fbdata%plam(jo) = profdata%rlam(jo) |
---|
659 | fbdata%pphi(jo) = profdata%rphi(jo) |
---|
660 | WRITE(fbdata%cdtyp(jo),'(I4)') profdata%ntyp(jo) |
---|
661 | fbdata%ivqc(jo,:) = profdata%ivqc(jo,:) |
---|
662 | fbdata%ivqcf(:,jo,:) = profdata%ivqcf(:,jo,:) |
---|
663 | IF ( profdata%nqc(jo) > 10 ) THEN |
---|
664 | fbdata%ioqc(jo) = 4 |
---|
665 | fbdata%ioqcf(1,jo) = profdata%nqcf(1,jo) |
---|
666 | fbdata%ioqcf(2,jo) = profdata%nqc(jo) - 10 |
---|
667 | ELSE |
---|
668 | fbdata%ioqc(jo) = profdata%nqc(jo) |
---|
669 | fbdata%ioqcf(:,jo) = profdata%nqcf(:,jo) |
---|
670 | ENDIF |
---|
671 | fbdata%ipqc(jo) = profdata%ipqc(jo) |
---|
672 | fbdata%ipqcf(:,jo) = profdata%ipqcf(:,jo) |
---|
673 | fbdata%itqc(jo) = profdata%itqc(jo) |
---|
674 | fbdata%itqcf(:,jo) = profdata%itqcf(:,jo) |
---|
675 | fbdata%cdwmo(jo) = profdata%cwmo(jo) |
---|
676 | fbdata%kindex(jo) = profdata%npfil(jo) |
---|
677 | DO jvar = 1, profdata%nvar |
---|
678 | IF (ln_grid_global) THEN |
---|
679 | fbdata%iobsi(jo,jvar) = profdata%mi(jo,jvar) |
---|
680 | fbdata%iobsj(jo,jvar) = profdata%mj(jo,jvar) |
---|
681 | ELSE |
---|
682 | fbdata%iobsi(jo,jvar) = mig(profdata%mi(jo,jvar)) |
---|
683 | fbdata%iobsj(jo,jvar) = mjg(profdata%mj(jo,jvar)) |
---|
684 | ENDIF |
---|
685 | ENDDO |
---|
686 | CALL greg2jul( 0, & |
---|
687 | & profdata%nmin(jo), & |
---|
688 | & profdata%nhou(jo), & |
---|
689 | & profdata%nday(jo), & |
---|
690 | & profdata%nmon(jo), & |
---|
691 | & profdata%nyea(jo), & |
---|
692 | & fbdata%ptim(jo), & |
---|
693 | & krefdate = 19500101 ) |
---|
694 | ! Reform the profiles arrays for output |
---|
695 | DO jvar = 1, 2 |
---|
696 | DO jk = profdata%npvsta(jo,jvar), profdata%npvend(jo,jvar) |
---|
697 | ik = profdata%var(jvar)%nvlidx(jk) |
---|
698 | IF ( jvar == 1 ) THEN |
---|
699 | fbdata%padd(ik,jo,1,jvar) = zu(jk) |
---|
700 | ELSE |
---|
701 | fbdata%padd(ik,jo,1,jvar) = zv(jk) |
---|
702 | ENDIF |
---|
703 | fbdata%padd(ik,jo,2,jvar) = profdata%var(jvar)%vmod(jk) |
---|
704 | fbdata%pob(ik,jo,jvar) = profdata%var(jvar)%vobs(jk) |
---|
705 | fbdata%pdep(ik,jo) = profdata%var(jvar)%vdep(jk) |
---|
706 | fbdata%idqc(ik,jo) = profdata%var(jvar)%idqc(jk) |
---|
707 | fbdata%idqcf(:,ik,jo) = profdata%var(jvar)%idqcf(:,jk) |
---|
708 | IF ( profdata%var(jvar)%nvqc(jk) > 10 ) THEN |
---|
709 | fbdata%ivlqc(ik,jo,jvar) = 4 |
---|
710 | fbdata%ivlqcf(1,ik,jo,jvar) = profdata%var(jvar)%nvqcf(1,jk) |
---|
711 | fbdata%ivlqcf(2,ik,jo,jvar) = profdata%var(jvar)%nvqc(jk) - 10 |
---|
712 | ELSE |
---|
713 | fbdata%ivlqc(ik,jo,jvar) = profdata%var(jvar)%nvqc(jk) |
---|
714 | fbdata%ivlqcf(:,ik,jo,jvar) = profdata%var(jvar)%nvqcf(:,jk) |
---|
715 | ENDIF |
---|
716 | fbdata%iobsk(ik,jo,jvar) = profdata%var(jvar)%mvk(jk) |
---|
717 | END DO |
---|
718 | END DO |
---|
719 | END DO |
---|
720 | |
---|
721 | ! Write the obfbdata structure |
---|
722 | CALL write_obfbdata( cfname, fbdata ) |
---|
723 | |
---|
724 | CALL dealloc_obfbdata( fbdata ) |
---|
725 | |
---|
726 | DEALLOCATE( & |
---|
727 | & zu, & |
---|
728 | & zv & |
---|
729 | & ) |
---|
730 | |
---|
731 | END SUBROUTINE obs_wri_vel |
---|
732 | |
---|
733 | END MODULE obs_write |
---|