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 | & obswriinfo |
---|
45 | |
---|
46 | TYPE obswriinfo |
---|
47 | INTEGER :: inum |
---|
48 | INTEGER, POINTER, DIMENSION(:) :: ipoint |
---|
49 | CHARACTER(len=ilenname), POINTER, DIMENSION(:) :: cdname |
---|
50 | CHARACTER(len=ilenlong), POINTER, DIMENSION(:,:) :: cdlong |
---|
51 | CHARACTER(len=ilenunit), POINTER, DIMENSION(:,:) :: cdunit |
---|
52 | END TYPE obswriinfo |
---|
53 | |
---|
54 | !!---------------------------------------------------------------------- |
---|
55 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
56 | !! $Id$ |
---|
57 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
58 | !!---------------------------------------------------------------------- |
---|
59 | |
---|
60 | CONTAINS |
---|
61 | |
---|
62 | SUBROUTINE obs_wri_p3d( cprefix, profdata, padd, pext ) |
---|
63 | !!----------------------------------------------------------------------- |
---|
64 | !! |
---|
65 | !! *** ROUTINE obs_wri_p3d *** |
---|
66 | !! |
---|
67 | !! ** Purpose : Write temperature and salinity (profile) observation |
---|
68 | !! related diagnostics |
---|
69 | !! |
---|
70 | !! ** Method : NetCDF |
---|
71 | !! |
---|
72 | !! ** Action : |
---|
73 | !! |
---|
74 | !! History : |
---|
75 | !! ! 06-04 (A. Vidard) Original |
---|
76 | !! ! 06-04 (A. Vidard) Reformatted |
---|
77 | !! ! 06-10 (A. Weaver) Cleanup |
---|
78 | !! ! 07-01 (K. Mogensen) Use profile data types |
---|
79 | !! ! 07-03 (K. Mogensen) General handling of profiles |
---|
80 | !! ! 09-01 (K. Mogensen) New feedback format |
---|
81 | !!----------------------------------------------------------------------- |
---|
82 | |
---|
83 | !! * Modules used |
---|
84 | |
---|
85 | !! * Arguments |
---|
86 | CHARACTER(LEN=*), INTENT(IN) :: cprefix ! Prefix for output files |
---|
87 | TYPE(obs_prof), INTENT(INOUT) :: profdata ! Full set of profile data |
---|
88 | TYPE(obswriinfo), OPTIONAL :: padd ! Additional info for each variable |
---|
89 | TYPE(obswriinfo), OPTIONAL :: pext ! Extra info |
---|
90 | |
---|
91 | !! * Local declarations |
---|
92 | TYPE(obfbdata) :: fbdata |
---|
93 | CHARACTER(LEN=40) :: cfname |
---|
94 | INTEGER :: ilevel |
---|
95 | INTEGER :: jvar |
---|
96 | INTEGER :: jo |
---|
97 | INTEGER :: jk |
---|
98 | INTEGER :: ik |
---|
99 | INTEGER :: ja |
---|
100 | INTEGER :: je |
---|
101 | REAL(wp) :: zpres |
---|
102 | INTEGER :: nadd |
---|
103 | INTEGER :: next |
---|
104 | |
---|
105 | IF ( PRESENT( padd ) ) THEN |
---|
106 | nadd = padd%inum |
---|
107 | ELSE |
---|
108 | nadd = 0 |
---|
109 | ENDIF |
---|
110 | |
---|
111 | IF ( PRESENT( pext ) ) THEN |
---|
112 | next = pext%inum |
---|
113 | ELSE |
---|
114 | next = 0 |
---|
115 | ENDIF |
---|
116 | |
---|
117 | CALL init_obfbdata( fbdata ) |
---|
118 | |
---|
119 | ! Find maximum level |
---|
120 | ilevel = 0 |
---|
121 | DO jvar = 1, 2 |
---|
122 | ilevel = MAX( ilevel, MAXVAL( profdata%var(jvar)%nvlidx(:) ) ) |
---|
123 | END DO |
---|
124 | CALL alloc_obfbdata( fbdata, 2, profdata%nprof, ilevel, & |
---|
125 | & 1 + nadd, 1 + next, .TRUE. ) |
---|
126 | |
---|
127 | fbdata%cname(1) = 'POTM' |
---|
128 | fbdata%cname(2) = 'PSAL' |
---|
129 | fbdata%coblong(1) = 'Potential temperature' |
---|
130 | fbdata%coblong(2) = 'Practical salinity' |
---|
131 | fbdata%cobunit(1) = 'Degrees centigrade' |
---|
132 | fbdata%cobunit(2) = 'PSU' |
---|
133 | fbdata%cextname(1) = 'TEMP' |
---|
134 | fbdata%cextlong(1) = 'Insitu temperature' |
---|
135 | fbdata%cextunit(1) = 'Degrees centigrade' |
---|
136 | DO je = 1, next |
---|
137 | fbdata%cextname(1+je) = pext%cdname(je) |
---|
138 | fbdata%cextlong(1+je) = pext%cdlong(je,1) |
---|
139 | fbdata%cextunit(1+je) = pext%cdunit(je,1) |
---|
140 | END DO |
---|
141 | fbdata%caddname(1) = 'Hx' |
---|
142 | fbdata%caddlong(1,1) = 'Model interpolated potential temperature' |
---|
143 | fbdata%caddlong(1,2) = 'Model interpolated practical salinity' |
---|
144 | fbdata%caddunit(1,1) = 'Degrees centigrade' |
---|
145 | fbdata%caddunit(1,2) = 'PSU' |
---|
146 | fbdata%cgrid(:) = 'T' |
---|
147 | DO ja = 1, nadd |
---|
148 | fbdata%caddname(1+ja) = padd%cdname(ja) |
---|
149 | DO jvar = 1, 2 |
---|
150 | fbdata%caddlong(1+ja,jvar) = padd%cdlong(ja,jvar) |
---|
151 | fbdata%caddunit(1+ja,jvar) = padd%cdunit(ja,jvar) |
---|
152 | END DO |
---|
153 | END DO |
---|
154 | |
---|
155 | WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc |
---|
156 | |
---|
157 | IF(lwp) THEN |
---|
158 | WRITE(numout,*) |
---|
159 | WRITE(numout,*)'obs_wri_p3d :' |
---|
160 | WRITE(numout,*)'~~~~~~~~~~~~~' |
---|
161 | WRITE(numout,*)'Writing profile feedback file : ',TRIM(cfname) |
---|
162 | ENDIF |
---|
163 | |
---|
164 | ! Transform obs_prof data structure into obfbdata structure |
---|
165 | fbdata%cdjuldref = '19500101000000' |
---|
166 | DO jo = 1, profdata%nprof |
---|
167 | fbdata%plam(jo) = profdata%rlam(jo) |
---|
168 | fbdata%pphi(jo) = profdata%rphi(jo) |
---|
169 | WRITE(fbdata%cdtyp(jo),'(I4)') profdata%ntyp(jo) |
---|
170 | fbdata%ivqc(jo,:) = profdata%ivqc(jo,:) |
---|
171 | fbdata%ivqcf(:,jo,:) = profdata%ivqcf(:,jo,:) |
---|
172 | IF ( profdata%nqc(jo) > 10 ) THEN |
---|
173 | fbdata%ioqc(jo) = 4 |
---|
174 | fbdata%ioqcf(1,jo) = profdata%nqcf(1,jo) |
---|
175 | fbdata%ioqcf(2,jo) = profdata%nqc(jo) - 10 |
---|
176 | ELSE |
---|
177 | fbdata%ioqc(jo) = profdata%nqc(jo) |
---|
178 | fbdata%ioqcf(:,jo) = profdata%nqcf(:,jo) |
---|
179 | ENDIF |
---|
180 | fbdata%ipqc(jo) = profdata%ipqc(jo) |
---|
181 | fbdata%ipqcf(:,jo) = profdata%ipqcf(:,jo) |
---|
182 | fbdata%itqc(jo) = profdata%itqc(jo) |
---|
183 | fbdata%itqcf(:,jo) = profdata%itqcf(:,jo) |
---|
184 | fbdata%cdwmo(jo) = profdata%cwmo(jo) |
---|
185 | fbdata%kindex(jo) = profdata%npfil(jo) |
---|
186 | DO jvar = 1, profdata%nvar |
---|
187 | IF (ln_grid_global) THEN |
---|
188 | fbdata%iobsi(jo,jvar) = profdata%mi(jo,jvar) |
---|
189 | fbdata%iobsj(jo,jvar) = profdata%mj(jo,jvar) |
---|
190 | ELSE |
---|
191 | fbdata%iobsi(jo,jvar) = mig(profdata%mi(jo,jvar)) |
---|
192 | fbdata%iobsj(jo,jvar) = mjg(profdata%mj(jo,jvar)) |
---|
193 | ENDIF |
---|
194 | END DO |
---|
195 | CALL greg2jul( 0, & |
---|
196 | & profdata%nmin(jo), & |
---|
197 | & profdata%nhou(jo), & |
---|
198 | & profdata%nday(jo), & |
---|
199 | & profdata%nmon(jo), & |
---|
200 | & profdata%nyea(jo), & |
---|
201 | & fbdata%ptim(jo), & |
---|
202 | & krefdate = 19500101 ) |
---|
203 | ! Reform the profiles arrays for output |
---|
204 | DO jvar = 1, 2 |
---|
205 | DO jk = profdata%npvsta(jo,jvar), profdata%npvend(jo,jvar) |
---|
206 | ik = profdata%var(jvar)%nvlidx(jk) |
---|
207 | fbdata%padd(ik,jo,1,jvar) = profdata%var(jvar)%vmod(jk) |
---|
208 | fbdata%pob(ik,jo,jvar) = profdata%var(jvar)%vobs(jk) |
---|
209 | fbdata%pdep(ik,jo) = profdata%var(jvar)%vdep(jk) |
---|
210 | fbdata%idqc(ik,jo) = profdata%var(jvar)%idqc(jk) |
---|
211 | fbdata%idqcf(:,ik,jo) = profdata%var(jvar)%idqcf(:,jk) |
---|
212 | IF ( profdata%var(jvar)%nvqc(jk) > 10 ) THEN |
---|
213 | fbdata%ivlqc(ik,jo,jvar) = 4 |
---|
214 | fbdata%ivlqcf(1,ik,jo,jvar) = profdata%var(jvar)%nvqcf(1,jk) |
---|
215 | fbdata%ivlqcf(2,ik,jo,jvar) = profdata%var(jvar)%nvqc(jk) - 10 |
---|
216 | ELSE |
---|
217 | fbdata%ivlqc(ik,jo,jvar) = profdata%var(jvar)%nvqc(jk) |
---|
218 | fbdata%ivlqcf(:,ik,jo,jvar) = profdata%var(jvar)%nvqcf(:,jk) |
---|
219 | ENDIF |
---|
220 | fbdata%iobsk(ik,jo,jvar) = profdata%var(jvar)%mvk(jk) |
---|
221 | DO ja = 1, nadd |
---|
222 | fbdata%padd(ik,jo,1+ja,jvar) = & |
---|
223 | & profdata%var(jvar)%vext(jk,padd%ipoint(ja)) |
---|
224 | END DO |
---|
225 | DO je = 1, next |
---|
226 | fbdata%pext(ik,jo,1+je) = & |
---|
227 | & profdata%var(jvar)%vext(jk,pext%ipoint(je)) |
---|
228 | END DO |
---|
229 | IF ( jvar == 1 ) THEN |
---|
230 | fbdata%pext(ik,jo,1) = profdata%var(jvar)%vext(jk,1) |
---|
231 | ENDIF |
---|
232 | END DO |
---|
233 | END DO |
---|
234 | END DO |
---|
235 | |
---|
236 | ! Convert insitu temperature to potential temperature using the model |
---|
237 | ! salinity if no potential temperature |
---|
238 | DO jo = 1, fbdata%nobs |
---|
239 | IF ( fbdata%pphi(jo) < 9999.0 ) THEN |
---|
240 | DO jk = 1, fbdata%nlev |
---|
241 | IF ( ( fbdata%pob(jk,jo,1) >= 9999.0 ) .AND. & |
---|
242 | & ( fbdata%pdep(jk,jo) < 9999.0 ) .AND. & |
---|
243 | & ( fbdata%padd(jk,jo,1,2) < 9999.0 ) .AND. & |
---|
244 | & ( fbdata%pext(jk,jo,1) < 9999.0 ) ) THEN |
---|
245 | zpres = dep_to_p( REAL(fbdata%pdep(jk,jo),wp), & |
---|
246 | & REAL(fbdata%pphi(jo),wp) ) |
---|
247 | fbdata%pob(jk,jo,1) = potemp( & |
---|
248 | & REAL(fbdata%padd(jk,jo,1,2), wp), & |
---|
249 | & REAL(fbdata%pext(jk,jo,1), wp), & |
---|
250 | & zpres, 0.0_wp ) |
---|
251 | ENDIF |
---|
252 | END DO |
---|
253 | ENDIF |
---|
254 | END DO |
---|
255 | |
---|
256 | ! Write the obfbdata structure |
---|
257 | CALL write_obfbdata( cfname, fbdata ) |
---|
258 | |
---|
259 | CALL dealloc_obfbdata( fbdata ) |
---|
260 | |
---|
261 | END SUBROUTINE obs_wri_p3d |
---|
262 | |
---|
263 | SUBROUTINE obs_wri_sla( cprefix, sladata, padd, pext ) |
---|
264 | !!----------------------------------------------------------------------- |
---|
265 | !! |
---|
266 | !! *** ROUTINE obs_wri_sla *** |
---|
267 | !! |
---|
268 | !! ** Purpose : Write SLA observation diagnostics |
---|
269 | !! related |
---|
270 | !! |
---|
271 | !! ** Method : NetCDF |
---|
272 | !! |
---|
273 | !! ** Action : |
---|
274 | !! |
---|
275 | !! ! 07-03 (K. Mogensen) Original |
---|
276 | !! ! 09-01 (K. Mogensen) New feedback format. |
---|
277 | !!----------------------------------------------------------------------- |
---|
278 | |
---|
279 | !! * Modules used |
---|
280 | IMPLICIT NONE |
---|
281 | |
---|
282 | !! * Arguments |
---|
283 | CHARACTER(LEN=*), INTENT(IN) :: cprefix ! Prefix for output files |
---|
284 | TYPE(obs_surf), INTENT(INOUT) :: sladata ! Full set of SLAa |
---|
285 | TYPE(obswriinfo), OPTIONAL :: padd ! Additional info for each variable |
---|
286 | TYPE(obswriinfo), OPTIONAL :: pext ! Extra info |
---|
287 | |
---|
288 | !! * Local declarations |
---|
289 | TYPE(obfbdata) :: fbdata |
---|
290 | CHARACTER(LEN=40) :: cfname ! netCDF filename |
---|
291 | CHARACTER(LEN=12), PARAMETER :: cpname = 'obs_wri_sla' |
---|
292 | INTEGER :: jo |
---|
293 | INTEGER :: ja |
---|
294 | INTEGER :: je |
---|
295 | INTEGER :: nadd |
---|
296 | INTEGER :: next |
---|
297 | |
---|
298 | IF ( PRESENT( padd ) ) THEN |
---|
299 | nadd = padd%inum |
---|
300 | ELSE |
---|
301 | nadd = 0 |
---|
302 | ENDIF |
---|
303 | |
---|
304 | IF ( PRESENT( pext ) ) THEN |
---|
305 | next = pext%inum |
---|
306 | ELSE |
---|
307 | next = 0 |
---|
308 | ENDIF |
---|
309 | |
---|
310 | CALL init_obfbdata( fbdata ) |
---|
311 | |
---|
312 | CALL alloc_obfbdata( fbdata, 1, sladata%nsurf, 1, & |
---|
313 | & 2 + nadd, 1 + next, .TRUE. ) |
---|
314 | |
---|
315 | fbdata%cname(1) = 'SLA' |
---|
316 | fbdata%coblong(1) = 'Sea level anomaly' |
---|
317 | fbdata%cobunit(1) = 'Metres' |
---|
318 | fbdata%cextname(1) = 'MDT' |
---|
319 | fbdata%cextlong(1) = 'Mean dynamic topography' |
---|
320 | fbdata%cextunit(1) = 'Metres' |
---|
321 | DO je = 1, next |
---|
322 | fbdata%cextname(1+je) = pext%cdname(je) |
---|
323 | fbdata%cextlong(1+je) = pext%cdlong(je,1) |
---|
324 | fbdata%cextunit(1+je) = pext%cdunit(je,1) |
---|
325 | END DO |
---|
326 | fbdata%caddname(1) = 'Hx' |
---|
327 | fbdata%caddlong(1,1) = 'Model interpolated SSH - MDT' |
---|
328 | fbdata%caddunit(1,1) = 'Metres' |
---|
329 | fbdata%caddname(2) = 'SSH' |
---|
330 | fbdata%caddlong(2,1) = 'Model Sea surface height' |
---|
331 | fbdata%caddunit(2,1) = 'Metres' |
---|
332 | fbdata%cgrid(1) = 'T' |
---|
333 | DO ja = 1, nadd |
---|
334 | fbdata%caddname(2+ja) = padd%cdname(ja) |
---|
335 | fbdata%caddlong(2+ja,1) = padd%cdlong(ja,1) |
---|
336 | fbdata%caddunit(2+ja,1) = padd%cdunit(ja,1) |
---|
337 | END DO |
---|
338 | |
---|
339 | WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc |
---|
340 | |
---|
341 | IF(lwp) THEN |
---|
342 | WRITE(numout,*) |
---|
343 | WRITE(numout,*)'obs_wri_sla :' |
---|
344 | WRITE(numout,*)'~~~~~~~~~~~~~' |
---|
345 | WRITE(numout,*)'Writing SLA feedback file : ',TRIM(cfname) |
---|
346 | ENDIF |
---|
347 | |
---|
348 | ! Transform obs_prof data structure into obfbdata structure |
---|
349 | fbdata%cdjuldref = '19500101000000' |
---|
350 | DO jo = 1, sladata%nsurf |
---|
351 | fbdata%plam(jo) = sladata%rlam(jo) |
---|
352 | fbdata%pphi(jo) = sladata%rphi(jo) |
---|
353 | WRITE(fbdata%cdtyp(jo),'(I4)') sladata%ntyp(jo) |
---|
354 | fbdata%ivqc(jo,:) = 0 |
---|
355 | fbdata%ivqcf(:,jo,:) = 0 |
---|
356 | IF ( sladata%nqc(jo) > 10 ) THEN |
---|
357 | fbdata%ioqc(jo) = 4 |
---|
358 | fbdata%ioqcf(1,jo) = 0 |
---|
359 | fbdata%ioqcf(2,jo) = sladata%nqc(jo) - 10 |
---|
360 | ELSE |
---|
361 | fbdata%ioqc(jo) = sladata%nqc(jo) |
---|
362 | fbdata%ioqcf(:,jo) = 0 |
---|
363 | ENDIF |
---|
364 | fbdata%ipqc(jo) = 0 |
---|
365 | fbdata%ipqcf(:,jo) = 0 |
---|
366 | fbdata%itqc(jo) = 0 |
---|
367 | fbdata%itqcf(:,jo) = 0 |
---|
368 | fbdata%cdwmo(jo) = sladata%cwmo(jo) |
---|
369 | fbdata%kindex(jo) = sladata%nsfil(jo) |
---|
370 | IF (ln_grid_global) THEN |
---|
371 | fbdata%iobsi(jo,1) = sladata%mi(jo) |
---|
372 | fbdata%iobsj(jo,1) = sladata%mj(jo) |
---|
373 | ELSE |
---|
374 | fbdata%iobsi(jo,1) = mig(sladata%mi(jo)) |
---|
375 | fbdata%iobsj(jo,1) = mjg(sladata%mj(jo)) |
---|
376 | ENDIF |
---|
377 | CALL greg2jul( 0, & |
---|
378 | & sladata%nmin(jo), & |
---|
379 | & sladata%nhou(jo), & |
---|
380 | & sladata%nday(jo), & |
---|
381 | & sladata%nmon(jo), & |
---|
382 | & sladata%nyea(jo), & |
---|
383 | & fbdata%ptim(jo), & |
---|
384 | & krefdate = 19500101 ) |
---|
385 | fbdata%padd(1,jo,1,1) = sladata%rmod(jo,1) |
---|
386 | fbdata%padd(1,jo,2,1) = sladata%rext(jo,1) |
---|
387 | fbdata%pob(1,jo,1) = sladata%robs(jo,1) |
---|
388 | fbdata%pdep(1,jo) = 0.0 |
---|
389 | fbdata%idqc(1,jo) = 0 |
---|
390 | fbdata%idqcf(:,1,jo) = 0 |
---|
391 | IF ( sladata%nqc(jo) > 10 ) THEN |
---|
392 | fbdata%ivqc(jo,1) = 4 |
---|
393 | fbdata%ivlqc(1,jo,1) = 4 |
---|
394 | fbdata%ivlqcf(1,1,jo,1) = 0 |
---|
395 | fbdata%ivlqcf(2,1,jo,1) = sladata%nqc(jo) - 10 |
---|
396 | ELSE |
---|
397 | fbdata%ivqc(jo,1) = sladata%nqc(jo) |
---|
398 | fbdata%ivlqc(1,jo,1) = sladata%nqc(jo) |
---|
399 | fbdata%ivlqcf(:,1,jo,1) = 0 |
---|
400 | ENDIF |
---|
401 | fbdata%iobsk(1,jo,1) = 0 |
---|
402 | fbdata%pext(1,jo,1) = sladata%rext(jo,2) |
---|
403 | DO ja = 1, nadd |
---|
404 | fbdata%padd(1,jo,2+ja,1) = & |
---|
405 | & sladata%rext(jo,padd%ipoint(ja)) |
---|
406 | END DO |
---|
407 | DO je = 1, next |
---|
408 | fbdata%pext(1,jo,1+je) = & |
---|
409 | & sladata%rext(jo,pext%ipoint(je)) |
---|
410 | END DO |
---|
411 | END DO |
---|
412 | |
---|
413 | ! Write the obfbdata structure |
---|
414 | CALL write_obfbdata( cfname, fbdata ) |
---|
415 | |
---|
416 | CALL dealloc_obfbdata( fbdata ) |
---|
417 | |
---|
418 | END SUBROUTINE obs_wri_sla |
---|
419 | |
---|
420 | SUBROUTINE obs_wri_sst( cprefix, sstdata, padd, pext ) |
---|
421 | !!----------------------------------------------------------------------- |
---|
422 | !! |
---|
423 | !! *** ROUTINE obs_wri_sst *** |
---|
424 | !! |
---|
425 | !! ** Purpose : Write SST observation diagnostics |
---|
426 | !! related |
---|
427 | !! |
---|
428 | !! ** Method : NetCDF |
---|
429 | !! |
---|
430 | !! ** Action : |
---|
431 | !! |
---|
432 | !! ! 07-07 (S. Ricci) Original |
---|
433 | !! ! 09-01 (K. Mogensen) New feedback format. |
---|
434 | !!----------------------------------------------------------------------- |
---|
435 | |
---|
436 | !! * Modules used |
---|
437 | IMPLICIT NONE |
---|
438 | |
---|
439 | !! * Arguments |
---|
440 | CHARACTER(LEN=*), INTENT(IN) :: cprefix ! Prefix for output files |
---|
441 | TYPE(obs_surf), INTENT(INOUT) :: sstdata ! Full set of SST |
---|
442 | TYPE(obswriinfo), OPTIONAL :: padd ! Additional info for each variable |
---|
443 | TYPE(obswriinfo), OPTIONAL :: pext ! Extra info |
---|
444 | |
---|
445 | !! * Local declarations |
---|
446 | TYPE(obfbdata) :: fbdata |
---|
447 | CHARACTER(LEN=40) :: cfname ! netCDF filename |
---|
448 | CHARACTER(LEN=12), PARAMETER :: cpname = 'obs_wri_sst' |
---|
449 | INTEGER :: jo |
---|
450 | INTEGER :: ja |
---|
451 | INTEGER :: je |
---|
452 | INTEGER :: nadd |
---|
453 | INTEGER :: next |
---|
454 | |
---|
455 | IF ( PRESENT( padd ) ) THEN |
---|
456 | nadd = padd%inum |
---|
457 | ELSE |
---|
458 | nadd = 0 |
---|
459 | ENDIF |
---|
460 | |
---|
461 | IF ( PRESENT( pext ) ) THEN |
---|
462 | next = pext%inum |
---|
463 | ELSE |
---|
464 | next = 0 |
---|
465 | ENDIF |
---|
466 | |
---|
467 | CALL init_obfbdata( fbdata ) |
---|
468 | |
---|
469 | CALL alloc_obfbdata( fbdata, 1, sstdata%nsurf, 1, & |
---|
470 | & 1 + nadd, next, .TRUE. ) |
---|
471 | |
---|
472 | fbdata%cname(1) = 'SST' |
---|
473 | fbdata%coblong(1) = 'Sea surface temperature' |
---|
474 | fbdata%cobunit(1) = 'Degree centigrade' |
---|
475 | DO je = 1, next |
---|
476 | fbdata%cextname(je) = pext%cdname(je) |
---|
477 | fbdata%cextlong(je) = pext%cdlong(je,1) |
---|
478 | fbdata%cextunit(je) = pext%cdunit(je,1) |
---|
479 | END DO |
---|
480 | fbdata%caddname(1) = 'Hx' |
---|
481 | fbdata%caddlong(1,1) = 'Model interpolated SST' |
---|
482 | fbdata%caddunit(1,1) = 'Degree centigrade' |
---|
483 | fbdata%cgrid(1) = 'T' |
---|
484 | DO ja = 1, nadd |
---|
485 | fbdata%caddname(1+ja) = padd%cdname(ja) |
---|
486 | fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1) |
---|
487 | fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1) |
---|
488 | END DO |
---|
489 | |
---|
490 | WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc |
---|
491 | |
---|
492 | IF(lwp) THEN |
---|
493 | WRITE(numout,*) |
---|
494 | WRITE(numout,*)'obs_wri_sst :' |
---|
495 | WRITE(numout,*)'~~~~~~~~~~~~~' |
---|
496 | WRITE(numout,*)'Writing SST feedback file : ',TRIM(cfname) |
---|
497 | ENDIF |
---|
498 | |
---|
499 | ! Transform obs_prof data structure into obfbdata structure |
---|
500 | fbdata%cdjuldref = '19500101000000' |
---|
501 | DO jo = 1, sstdata%nsurf |
---|
502 | fbdata%plam(jo) = sstdata%rlam(jo) |
---|
503 | fbdata%pphi(jo) = sstdata%rphi(jo) |
---|
504 | WRITE(fbdata%cdtyp(jo),'(I4)') sstdata%ntyp(jo) |
---|
505 | fbdata%ivqc(jo,:) = 0 |
---|
506 | fbdata%ivqcf(:,jo,:) = 0 |
---|
507 | IF ( sstdata%nqc(jo) > 10 ) THEN |
---|
508 | fbdata%ioqc(jo) = 4 |
---|
509 | fbdata%ioqcf(1,jo) = 0 |
---|
510 | fbdata%ioqcf(2,jo) = sstdata%nqc(jo) - 10 |
---|
511 | ELSE |
---|
512 | fbdata%ioqc(jo) = MAX(sstdata%nqc(jo),1) |
---|
513 | fbdata%ioqcf(:,jo) = 0 |
---|
514 | ENDIF |
---|
515 | fbdata%ipqc(jo) = 0 |
---|
516 | fbdata%ipqcf(:,jo) = 0 |
---|
517 | fbdata%itqc(jo) = 0 |
---|
518 | fbdata%itqcf(:,jo) = 0 |
---|
519 | fbdata%cdwmo(jo) = '' |
---|
520 | fbdata%kindex(jo) = sstdata%nsfil(jo) |
---|
521 | IF (ln_grid_global) THEN |
---|
522 | fbdata%iobsi(jo,1) = sstdata%mi(jo) |
---|
523 | fbdata%iobsj(jo,1) = sstdata%mj(jo) |
---|
524 | ELSE |
---|
525 | fbdata%iobsi(jo,1) = mig(sstdata%mi(jo)) |
---|
526 | fbdata%iobsj(jo,1) = mjg(sstdata%mj(jo)) |
---|
527 | ENDIF |
---|
528 | CALL greg2jul( 0, & |
---|
529 | & sstdata%nmin(jo), & |
---|
530 | & sstdata%nhou(jo), & |
---|
531 | & sstdata%nday(jo), & |
---|
532 | & sstdata%nmon(jo), & |
---|
533 | & sstdata%nyea(jo), & |
---|
534 | & fbdata%ptim(jo), & |
---|
535 | & krefdate = 19500101 ) |
---|
536 | fbdata%padd(1,jo,1,1) = sstdata%rmod(jo,1) |
---|
537 | fbdata%pob(1,jo,1) = sstdata%robs(jo,1) |
---|
538 | fbdata%pdep(1,jo) = 0.0 |
---|
539 | fbdata%idqc(1,jo) = 0 |
---|
540 | fbdata%idqcf(:,1,jo) = 0 |
---|
541 | IF ( sstdata%nqc(jo) > 10 ) THEN |
---|
542 | fbdata%ivqc(jo,1) = 4 |
---|
543 | fbdata%ivlqc(1,jo,1) = 4 |
---|
544 | fbdata%ivlqcf(1,1,jo,1) = 0 |
---|
545 | fbdata%ivlqcf(2,1,jo,1) = sstdata%nqc(jo) - 10 |
---|
546 | ELSE |
---|
547 | fbdata%ivqc(jo,1) = MAX(sstdata%nqc(jo),1) |
---|
548 | fbdata%ivlqc(1,jo,1) = MAX(sstdata%nqc(jo),1) |
---|
549 | fbdata%ivlqcf(:,1,jo,1) = 0 |
---|
550 | ENDIF |
---|
551 | fbdata%iobsk(1,jo,1) = 0 |
---|
552 | DO ja = 1, nadd |
---|
553 | fbdata%padd(1,jo,1+ja,1) = & |
---|
554 | & sstdata%rext(jo,padd%ipoint(ja)) |
---|
555 | END DO |
---|
556 | DO je = 1, next |
---|
557 | fbdata%pext(1,jo,je) = & |
---|
558 | & sstdata%rext(jo,pext%ipoint(je)) |
---|
559 | END DO |
---|
560 | |
---|
561 | END DO |
---|
562 | |
---|
563 | ! Write the obfbdata structure |
---|
564 | |
---|
565 | CALL write_obfbdata( cfname, fbdata ) |
---|
566 | |
---|
567 | CALL dealloc_obfbdata( fbdata ) |
---|
568 | |
---|
569 | END SUBROUTINE obs_wri_sst |
---|
570 | |
---|
571 | SUBROUTINE obs_wri_sss |
---|
572 | END SUBROUTINE obs_wri_sss |
---|
573 | |
---|
574 | SUBROUTINE obs_wri_seaice( cprefix, seaicedata, padd, pext ) |
---|
575 | !!----------------------------------------------------------------------- |
---|
576 | !! |
---|
577 | !! *** ROUTINE obs_wri_seaice *** |
---|
578 | !! |
---|
579 | !! ** Purpose : Write sea ice observation diagnostics |
---|
580 | !! related |
---|
581 | !! |
---|
582 | !! ** Method : NetCDF |
---|
583 | !! |
---|
584 | !! ** Action : |
---|
585 | !! |
---|
586 | !! ! 07-07 (S. Ricci) Original |
---|
587 | !! ! 09-01 (K. Mogensen) New feedback format. |
---|
588 | !!----------------------------------------------------------------------- |
---|
589 | |
---|
590 | !! * Modules used |
---|
591 | IMPLICIT NONE |
---|
592 | |
---|
593 | !! * Arguments |
---|
594 | CHARACTER(LEN=*), INTENT(IN) :: cprefix ! Prefix for output files |
---|
595 | TYPE(obs_surf), INTENT(INOUT) :: seaicedata ! Full set of sea ice |
---|
596 | TYPE(obswriinfo), OPTIONAL :: padd ! Additional info for each variable |
---|
597 | TYPE(obswriinfo), OPTIONAL :: pext ! Extra info |
---|
598 | |
---|
599 | !! * Local declarations |
---|
600 | TYPE(obfbdata) :: fbdata |
---|
601 | CHARACTER(LEN=40) :: cfname ! netCDF filename |
---|
602 | CHARACTER(LEN=12), PARAMETER :: cpname = 'obs_wri_seaice' |
---|
603 | INTEGER :: jo |
---|
604 | INTEGER :: ja |
---|
605 | INTEGER :: je |
---|
606 | INTEGER :: nadd |
---|
607 | INTEGER :: next |
---|
608 | |
---|
609 | IF ( PRESENT( padd ) ) THEN |
---|
610 | nadd = padd%inum |
---|
611 | ELSE |
---|
612 | nadd = 0 |
---|
613 | ENDIF |
---|
614 | |
---|
615 | IF ( PRESENT( pext ) ) THEN |
---|
616 | next = pext%inum |
---|
617 | ELSE |
---|
618 | next = 0 |
---|
619 | ENDIF |
---|
620 | |
---|
621 | CALL init_obfbdata( fbdata ) |
---|
622 | |
---|
623 | CALL alloc_obfbdata( fbdata, 1, seaicedata%nsurf, 1, 1, 0, .TRUE. ) |
---|
624 | |
---|
625 | fbdata%cname(1) = 'SEAICE' |
---|
626 | fbdata%coblong(1) = 'Sea ice' |
---|
627 | fbdata%cobunit(1) = 'Fraction' |
---|
628 | DO je = 1, next |
---|
629 | fbdata%cextname(je) = pext%cdname(je) |
---|
630 | fbdata%cextlong(je) = pext%cdlong(je,1) |
---|
631 | fbdata%cextunit(je) = pext%cdunit(je,1) |
---|
632 | END DO |
---|
633 | fbdata%caddname(1) = 'Hx' |
---|
634 | fbdata%caddlong(1,1) = 'Model interpolated ICE' |
---|
635 | fbdata%caddunit(1,1) = 'Fraction' |
---|
636 | fbdata%cgrid(1) = 'T' |
---|
637 | DO ja = 1, nadd |
---|
638 | fbdata%caddname(1+ja) = padd%cdname(ja) |
---|
639 | fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1) |
---|
640 | fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1) |
---|
641 | END DO |
---|
642 | |
---|
643 | WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc |
---|
644 | |
---|
645 | IF(lwp) THEN |
---|
646 | WRITE(numout,*) |
---|
647 | WRITE(numout,*)'obs_wri_seaice :' |
---|
648 | WRITE(numout,*)'~~~~~~~~~~~~~~~~' |
---|
649 | WRITE(numout,*)'Writing SEAICE feedback file : ',TRIM(cfname) |
---|
650 | ENDIF |
---|
651 | |
---|
652 | ! Transform obs_prof data structure into obfbdata structure |
---|
653 | fbdata%cdjuldref = '19500101000000' |
---|
654 | DO jo = 1, seaicedata%nsurf |
---|
655 | fbdata%plam(jo) = seaicedata%rlam(jo) |
---|
656 | fbdata%pphi(jo) = seaicedata%rphi(jo) |
---|
657 | WRITE(fbdata%cdtyp(jo),'(I4)') seaicedata%ntyp(jo) |
---|
658 | fbdata%ivqc(jo,:) = 0 |
---|
659 | fbdata%ivqcf(:,jo,:) = 0 |
---|
660 | IF ( seaicedata%nqc(jo) > 10 ) THEN |
---|
661 | fbdata%ioqc(jo) = 4 |
---|
662 | fbdata%ioqcf(1,jo) = 0 |
---|
663 | fbdata%ioqcf(2,jo) = seaicedata%nqc(jo) - 10 |
---|
664 | ELSE |
---|
665 | fbdata%ioqc(jo) = MAX(seaicedata%nqc(jo),1) |
---|
666 | fbdata%ioqcf(:,jo) = 0 |
---|
667 | ENDIF |
---|
668 | fbdata%ipqc(jo) = 0 |
---|
669 | fbdata%ipqcf(:,jo) = 0 |
---|
670 | fbdata%itqc(jo) = 0 |
---|
671 | fbdata%itqcf(:,jo) = 0 |
---|
672 | fbdata%cdwmo(jo) = '' |
---|
673 | fbdata%kindex(jo) = seaicedata%nsfil(jo) |
---|
674 | IF (ln_grid_global) THEN |
---|
675 | fbdata%iobsi(jo,1) = seaicedata%mi(jo) |
---|
676 | fbdata%iobsj(jo,1) = seaicedata%mj(jo) |
---|
677 | ELSE |
---|
678 | fbdata%iobsi(jo,1) = mig(seaicedata%mi(jo)) |
---|
679 | fbdata%iobsj(jo,1) = mjg(seaicedata%mj(jo)) |
---|
680 | ENDIF |
---|
681 | CALL greg2jul( 0, & |
---|
682 | & seaicedata%nmin(jo), & |
---|
683 | & seaicedata%nhou(jo), & |
---|
684 | & seaicedata%nday(jo), & |
---|
685 | & seaicedata%nmon(jo), & |
---|
686 | & seaicedata%nyea(jo), & |
---|
687 | & fbdata%ptim(jo), & |
---|
688 | & krefdate = 19500101 ) |
---|
689 | fbdata%padd(1,jo,1,1) = seaicedata%rmod(jo,1) |
---|
690 | fbdata%pob(1,jo,1) = seaicedata%robs(jo,1) |
---|
691 | fbdata%pdep(1,jo) = 0.0 |
---|
692 | fbdata%idqc(1,jo) = 0 |
---|
693 | fbdata%idqcf(:,1,jo) = 0 |
---|
694 | IF ( seaicedata%nqc(jo) > 10 ) THEN |
---|
695 | fbdata%ivlqc(1,jo,1) = 4 |
---|
696 | fbdata%ivlqcf(1,1,jo,1) = 0 |
---|
697 | fbdata%ivlqcf(2,1,jo,1) = seaicedata%nqc(jo) - 10 |
---|
698 | ELSE |
---|
699 | fbdata%ivlqc(1,jo,1) = MAX(seaicedata%nqc(jo),1) |
---|
700 | fbdata%ivlqcf(:,1,jo,1) = 0 |
---|
701 | ENDIF |
---|
702 | fbdata%iobsk(1,jo,1) = 0 |
---|
703 | DO ja = 1, nadd |
---|
704 | fbdata%padd(1,jo,1+ja,1) = & |
---|
705 | & seaicedata%rext(jo,padd%ipoint(ja)) |
---|
706 | END DO |
---|
707 | DO je = 1, next |
---|
708 | fbdata%pext(1,jo,je) = & |
---|
709 | & seaicedata%rext(jo,pext%ipoint(je)) |
---|
710 | END DO |
---|
711 | |
---|
712 | END DO |
---|
713 | |
---|
714 | ! Write the obfbdata structure |
---|
715 | CALL write_obfbdata( cfname, fbdata ) |
---|
716 | |
---|
717 | CALL dealloc_obfbdata( fbdata ) |
---|
718 | |
---|
719 | END SUBROUTINE obs_wri_seaice |
---|
720 | |
---|
721 | SUBROUTINE obs_wri_vel( cprefix, profdata, k2dint, padd, pext ) |
---|
722 | !!----------------------------------------------------------------------- |
---|
723 | !! |
---|
724 | !! *** ROUTINE obs_wri_p3d *** |
---|
725 | !! |
---|
726 | !! ** Purpose : Write current (profile) observation |
---|
727 | !! related diagnostics |
---|
728 | !! |
---|
729 | !! ** Method : NetCDF |
---|
730 | !! |
---|
731 | !! ** Action : |
---|
732 | !! |
---|
733 | !! History : |
---|
734 | !! ! 09-01 (K. Mogensen) New feedback format routine |
---|
735 | !!----------------------------------------------------------------------- |
---|
736 | |
---|
737 | !! * Modules used |
---|
738 | |
---|
739 | !! * Arguments |
---|
740 | CHARACTER(LEN=*), INTENT(IN) :: cprefix ! Prefix for output files |
---|
741 | TYPE(obs_prof), INTENT(INOUT) :: profdata ! Full set of profile data |
---|
742 | INTEGER, INTENT(IN) :: k2dint ! Horizontal interpolation method |
---|
743 | TYPE(obswriinfo), OPTIONAL :: padd ! Additional info for each variable |
---|
744 | TYPE(obswriinfo), OPTIONAL :: pext ! Extra info |
---|
745 | |
---|
746 | !! * Local declarations |
---|
747 | TYPE(obfbdata) :: fbdata |
---|
748 | CHARACTER(LEN=40) :: cfname |
---|
749 | INTEGER :: ilevel |
---|
750 | INTEGER :: jvar |
---|
751 | INTEGER :: jk |
---|
752 | INTEGER :: ik |
---|
753 | INTEGER :: jo |
---|
754 | INTEGER :: ja |
---|
755 | INTEGER :: je |
---|
756 | INTEGER :: nadd |
---|
757 | INTEGER :: next |
---|
758 | REAL(wp) :: zpres |
---|
759 | REAL(wp), DIMENSION(:), ALLOCATABLE :: & |
---|
760 | & zu, & |
---|
761 | & zv |
---|
762 | |
---|
763 | IF ( PRESENT( padd ) ) THEN |
---|
764 | nadd = padd%inum |
---|
765 | ELSE |
---|
766 | nadd = 0 |
---|
767 | ENDIF |
---|
768 | |
---|
769 | IF ( PRESENT( pext ) ) THEN |
---|
770 | next = pext%inum |
---|
771 | ELSE |
---|
772 | next = 0 |
---|
773 | ENDIF |
---|
774 | |
---|
775 | CALL init_obfbdata( fbdata ) |
---|
776 | |
---|
777 | ! Find maximum level |
---|
778 | ilevel = 0 |
---|
779 | DO jvar = 1, 2 |
---|
780 | ilevel = MAX( ilevel, MAXVAL( profdata%var(jvar)%nvlidx(:) ) ) |
---|
781 | END DO |
---|
782 | CALL alloc_obfbdata( fbdata, 2, profdata%nprof, ilevel, 2, 0, .TRUE. ) |
---|
783 | |
---|
784 | fbdata%cname(1) = 'UVEL' |
---|
785 | fbdata%cname(2) = 'VVEL' |
---|
786 | fbdata%coblong(1) = 'Zonal velocity' |
---|
787 | fbdata%coblong(2) = 'Meridional velocity' |
---|
788 | fbdata%cobunit(1) = 'm/s' |
---|
789 | fbdata%cobunit(2) = 'm/s' |
---|
790 | DO je = 1, next |
---|
791 | fbdata%cextname(je) = pext%cdname(je) |
---|
792 | fbdata%cextlong(je) = pext%cdlong(je,1) |
---|
793 | fbdata%cextunit(je) = pext%cdunit(je,1) |
---|
794 | END DO |
---|
795 | fbdata%caddname(1) = 'Hx' |
---|
796 | fbdata%caddlong(1,1) = 'Model interpolated zonal velocity' |
---|
797 | fbdata%caddlong(1,2) = 'Model interpolated meridional velocity' |
---|
798 | fbdata%caddunit(1,1) = 'm/s' |
---|
799 | fbdata%caddunit(1,2) = 'm/s' |
---|
800 | fbdata%caddname(2) = 'HxG' |
---|
801 | fbdata%caddlong(2,1) = 'Model interpolated zonal velocity (model grid)' |
---|
802 | fbdata%caddlong(2,2) = 'Model interpolated meridional velocity (model grid)' |
---|
803 | fbdata%caddunit(2,1) = 'm/s' |
---|
804 | fbdata%caddunit(2,2) = 'm/s' |
---|
805 | fbdata%cgrid(1) = 'U' |
---|
806 | fbdata%cgrid(2) = 'V' |
---|
807 | DO ja = 1, nadd |
---|
808 | fbdata%caddname(2+ja) = padd%cdname(ja) |
---|
809 | fbdata%caddlong(2+ja,1) = padd%cdlong(ja,1) |
---|
810 | fbdata%caddunit(2+ja,1) = padd%cdunit(ja,1) |
---|
811 | END DO |
---|
812 | |
---|
813 | WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc |
---|
814 | |
---|
815 | IF(lwp) THEN |
---|
816 | WRITE(numout,*) |
---|
817 | WRITE(numout,*)'obs_wri_vel :' |
---|
818 | WRITE(numout,*)'~~~~~~~~~~~~~' |
---|
819 | WRITE(numout,*)'Writing velocuty feedback file : ',TRIM(cfname) |
---|
820 | ENDIF |
---|
821 | |
---|
822 | ALLOCATE( & |
---|
823 | & zu(profdata%nvprot(1)), & |
---|
824 | & zv(profdata%nvprot(2)) & |
---|
825 | & ) |
---|
826 | CALL obs_rotvel( profdata, k2dint, zu, zv ) |
---|
827 | |
---|
828 | ! Transform obs_prof data structure into obfbdata structure |
---|
829 | fbdata%cdjuldref = '19500101000000' |
---|
830 | DO jo = 1, profdata%nprof |
---|
831 | fbdata%plam(jo) = profdata%rlam(jo) |
---|
832 | fbdata%pphi(jo) = profdata%rphi(jo) |
---|
833 | WRITE(fbdata%cdtyp(jo),'(I4)') profdata%ntyp(jo) |
---|
834 | fbdata%ivqc(jo,:) = profdata%ivqc(jo,:) |
---|
835 | fbdata%ivqcf(:,jo,:) = profdata%ivqcf(:,jo,:) |
---|
836 | IF ( profdata%nqc(jo) > 10 ) THEN |
---|
837 | fbdata%ioqc(jo) = 4 |
---|
838 | fbdata%ioqcf(1,jo) = profdata%nqcf(1,jo) |
---|
839 | fbdata%ioqcf(2,jo) = profdata%nqc(jo) - 10 |
---|
840 | ELSE |
---|
841 | fbdata%ioqc(jo) = profdata%nqc(jo) |
---|
842 | fbdata%ioqcf(:,jo) = profdata%nqcf(:,jo) |
---|
843 | ENDIF |
---|
844 | fbdata%ipqc(jo) = profdata%ipqc(jo) |
---|
845 | fbdata%ipqcf(:,jo) = profdata%ipqcf(:,jo) |
---|
846 | fbdata%itqc(jo) = profdata%itqc(jo) |
---|
847 | fbdata%itqcf(:,jo) = profdata%itqcf(:,jo) |
---|
848 | fbdata%cdwmo(jo) = profdata%cwmo(jo) |
---|
849 | fbdata%kindex(jo) = profdata%npfil(jo) |
---|
850 | DO jvar = 1, profdata%nvar |
---|
851 | IF (ln_grid_global) THEN |
---|
852 | fbdata%iobsi(jo,jvar) = profdata%mi(jo,jvar) |
---|
853 | fbdata%iobsj(jo,jvar) = profdata%mj(jo,jvar) |
---|
854 | ELSE |
---|
855 | fbdata%iobsi(jo,jvar) = mig(profdata%mi(jo,jvar)) |
---|
856 | fbdata%iobsj(jo,jvar) = mjg(profdata%mj(jo,jvar)) |
---|
857 | ENDIF |
---|
858 | END DO |
---|
859 | CALL greg2jul( 0, & |
---|
860 | & profdata%nmin(jo), & |
---|
861 | & profdata%nhou(jo), & |
---|
862 | & profdata%nday(jo), & |
---|
863 | & profdata%nmon(jo), & |
---|
864 | & profdata%nyea(jo), & |
---|
865 | & fbdata%ptim(jo), & |
---|
866 | & krefdate = 19500101 ) |
---|
867 | ! Reform the profiles arrays for output |
---|
868 | DO jvar = 1, 2 |
---|
869 | DO jk = profdata%npvsta(jo,jvar), profdata%npvend(jo,jvar) |
---|
870 | ik = profdata%var(jvar)%nvlidx(jk) |
---|
871 | IF ( jvar == 1 ) THEN |
---|
872 | fbdata%padd(ik,jo,1,jvar) = zu(jk) |
---|
873 | ELSE |
---|
874 | fbdata%padd(ik,jo,1,jvar) = zv(jk) |
---|
875 | ENDIF |
---|
876 | fbdata%padd(ik,jo,2,jvar) = profdata%var(jvar)%vmod(jk) |
---|
877 | fbdata%pob(ik,jo,jvar) = profdata%var(jvar)%vobs(jk) |
---|
878 | fbdata%pdep(ik,jo) = profdata%var(jvar)%vdep(jk) |
---|
879 | fbdata%idqc(ik,jo) = profdata%var(jvar)%idqc(jk) |
---|
880 | fbdata%idqcf(:,ik,jo) = profdata%var(jvar)%idqcf(:,jk) |
---|
881 | IF ( profdata%var(jvar)%nvqc(jk) > 10 ) THEN |
---|
882 | fbdata%ivlqc(ik,jo,jvar) = 4 |
---|
883 | fbdata%ivlqcf(1,ik,jo,jvar) = profdata%var(jvar)%nvqcf(1,jk) |
---|
884 | fbdata%ivlqcf(2,ik,jo,jvar) = profdata%var(jvar)%nvqc(jk) - 10 |
---|
885 | ELSE |
---|
886 | fbdata%ivlqc(ik,jo,jvar) = profdata%var(jvar)%nvqc(jk) |
---|
887 | fbdata%ivlqcf(:,ik,jo,jvar) = profdata%var(jvar)%nvqcf(:,jk) |
---|
888 | ENDIF |
---|
889 | fbdata%iobsk(ik,jo,jvar) = profdata%var(jvar)%mvk(jk) |
---|
890 | DO ja = 1, nadd |
---|
891 | fbdata%padd(ik,jo,2+ja,jvar) = & |
---|
892 | & profdata%var(jvar)%vext(jk,padd%ipoint(ja)) |
---|
893 | END DO |
---|
894 | DO je = 1, next |
---|
895 | fbdata%pext(ik,jo,je) = & |
---|
896 | & profdata%var(jvar)%vext(jk,pext%ipoint(je)) |
---|
897 | END DO |
---|
898 | END DO |
---|
899 | END DO |
---|
900 | END DO |
---|
901 | |
---|
902 | ! Write the obfbdata structure |
---|
903 | CALL write_obfbdata( cfname, fbdata ) |
---|
904 | |
---|
905 | CALL dealloc_obfbdata( fbdata ) |
---|
906 | |
---|
907 | DEALLOCATE( & |
---|
908 | & zu, & |
---|
909 | & zv & |
---|
910 | & ) |
---|
911 | |
---|
912 | END SUBROUTINE obs_wri_vel |
---|
913 | |
---|
914 | END MODULE obs_write |
---|