1 subroutine da_write_obs_etkf(ob, iv, re)
3 !-------------------------------------------------------------------------
4 ! Purpose: Writes out components of iv=O-B structure.
5 !-------------------------------------------------------------------------
7 ! Arthur P. Mizzi (NCAR/MMM) February 2011 Modfied to output the extended ob.etkf file.
11 type (y_type), intent(in) :: ob ! Observation structure.
12 type (iv_type), intent(in) :: iv ! O-B structure.
13 type (y_type), intent(inout) :: re ! residual vector.
15 integer :: n, k, num_obs, ios
16 integer :: ounit ! Output unit
17 character(len=20) :: filename
18 character(len=20) :: apm_char
19 integer :: apm_index, apm_int
22 if (trace_use) call da_trace_entry("da_write_obs_etkf")
24 !-------------------------------------------------------------------------
26 !-------------------------------------------------------------------------
30 call da_get_unit(ounit)
33 write(unit=filename, fmt='(a,i4.4)') 'ob.etkf.', myproc
35 write(unit=filename, fmt='(a)') 'ob.etkf.0000'
38 open (unit=ounit,file=trim(filename),form='formatted',status='replace', &
41 call da_error(__FILE__,__LINE__, &
42 (/"Cannot open ETKF observation file"//filename/))
45 ! [0] Format for extended ob.etkf files (APM 02-10-2011)
47 1000 format(3f17.7,2x,a10,2x,2(f6.0,2x),4(f8.2,2x),i10)
49 ! [1] Transfer surface obs:
51 if (iv%info(synop)%nlocal > 0) then
53 do n = 1, iv%info(synop)%nlocal
54 if (iv%info(synop)%proc_domain(1,n)) num_obs = num_obs + 1
58 do n = 1, iv%info(synop)%nlocal
59 if (iv%info(synop)%proc_domain(1,n)) then
62 ! read(iv%info(synop)%platform(n)(4:6),*) apm_int
63 if ( iv%synop(n)%u%qc >= 0 .and. ob%synop(n)%u /= missing_r ) then
64 apm_index = apm_index + 1
65 write(ounit,1000) ob%synop(n)%u, iv%synop(n)%u%inv, iv%synop(n)%u%error, &
66 'WNU', apm_plc, -888.88, iv%info(synop)%lat(1,n), iv%info(synop)%lon(1,n), &
67 (ob%synop(n)%p)/100., -888.88, apm_index
69 if ( iv%synop(n)%v%qc >= 0 .and. ob%synop(n)%v /= missing_r ) then
70 apm_index = apm_index + 1
71 write(ounit,1000) ob%synop(n)%v, iv%synop(n)%v%inv, iv%synop(n)%v%error, &
72 'WNV', apm_plc, -888.88, iv%info(synop)%lat(1,n), iv%info(synop)%lon(1,n), &
73 (ob%synop(n)%p)/100., -888.88, apm_index
75 if ( iv%synop(n)%t%qc >= 0 .and. ob%synop(n)%t /= missing_r ) then
76 apm_index = apm_index + 1
77 write(ounit,1000) ob%synop(n)%t, iv%synop(n)%t%inv, iv%synop(n)%t%error, &
78 'TMP', apm_plc, -888.88, iv%info(synop)%lat(1,n), iv%info(synop)%lon(1,n), &
79 (ob%synop(n)%p)/100., -888.88, apm_index
81 if ( iv%synop(n)%p%qc >= 0 .and. ob%synop(n)%p /= missing_r ) then
82 apm_index = apm_index + 1
83 write(ounit,1000) ob%synop(n)%p, iv%synop(n)%p%inv, iv%synop(n)%p%error, &
84 'PRS', apm_plc, -888.88, iv%info(synop)%lat(1,n), iv%info(synop)%lon(1,n), &
85 (ob%synop(n)%p)/100., -888.88, apm_index
87 if ( iv%synop(n)%q%qc >= 0 .and. ob%synop(n)%q /= missing_r ) then
88 apm_index = apm_index + 1
89 write(ounit,1000) ob%synop(n)%q, iv%synop(n)%q%inv, iv%synop(n)%q%error, &
90 'QVP', apm_plc, -888.88, iv%info(synop)%lat(1,n), iv%info(synop)%lon(1,n), &
91 (ob%synop(n)%p)/100., -888.88, apm_index
98 ! [2] Transfer metar obs:
100 if (iv%info(metar)%nlocal > 0) then
102 do n = 1, iv%info(metar)%nlocal
103 if (iv%info(metar)%proc_domain(1,n)) num_obs = num_obs + 1
105 if (num_obs > 0) then
107 do n = 1, iv%info(metar)%nlocal
108 if (iv%info(metar)%proc_domain(1,n)) then
109 num_obs = num_obs + 1
111 ! read(iv%info(metar)%platform(n)(4:6),*) apm_int
112 if ( iv%metar(n)%u%qc >= 0 .and. ob%metar(n)%u /= missing_r ) then
113 apm_index = apm_index + 1
114 write(ounit,1000) ob%metar(n)%u, iv%metar(n)%u%inv, iv%metar(n)%u%error, &
115 'WNU', apm_plc, -888.88, iv%info(metar)%lat(1,n), iv%info(metar)%lon(1,n), &
116 (ob%metar(n)%p)/100., -888.88, apm_index
118 if ( iv%metar(n)%v%qc >= 0 .and. ob%metar(n)%v /= missing_r ) then
119 apm_index = apm_index + 1
120 write(ounit,1000) ob%metar(n)%v, iv%metar(n)%v%inv, iv%metar(n)%v%error, &
121 'WNV', apm_plc, -888.88, iv%info(metar)%lat(1,n), iv%info(metar)%lon(1,n), &
122 (ob%metar(n)%p)/100., -888.88, apm_index
124 if ( iv%metar(n)%t%qc >= 0 .and. ob%metar(n)%t /= missing_r ) then
125 apm_index = apm_index + 1
126 write(ounit,1000) ob%metar(n)%t, iv%metar(n)%t%inv, iv%metar(n)%t%error, &
127 'TMP', apm_plc, -888.88, iv%info(metar)%lat(1,n), iv%info(metar)%lon(1,n), &
128 (ob%metar(n)%p)/100., -888.88, apm_index
130 if ( iv%metar(n)%p%qc >= 0 .and. ob%metar(n)%p /= missing_r ) then
131 apm_index = apm_index + 1
132 write(ounit,1000) ob%metar(n)%p, iv%metar(n)%p%inv, iv%metar(n)%p%error, &
133 'PRS', apm_plc, -888.88, iv%info(metar)%lat(1,n), iv%info(metar)%lon(1,n), &
134 (ob%metar(n)%p)/100., -888.88, apm_index
136 if ( iv%metar(n)%q%qc >= 0 .and. ob%metar(n)%q /= missing_r ) then
137 apm_index = apm_index + 1
138 write(ounit,1000) ob%metar(n)%q, iv%metar(n)%q%inv, iv%metar(n)%q%error, &
139 'QVP', apm_plc, -888.88, iv%info(metar)%lat(1,n), iv%info(metar)%lon(1,n), &
140 (ob%metar(n)%p)/100., -888.88, apm_index
147 ! [3] Transfer ships obs:
149 if (iv%info(ships)%nlocal > 0) then
151 do n = 1, iv%info(ships)%nlocal
152 if (iv%info(ships)%proc_domain(1,n)) num_obs = num_obs + 1
154 if (num_obs > 0) then
156 do n = 1, iv%info(ships)%nlocal
157 if (iv%info(ships)%proc_domain(1,n)) then
158 num_obs = num_obs + 1
160 ! read(iv%info(ships)%platform(n)(4:6),*) apm_int
161 if ( iv%ships(n)%u%qc >= 0 .and. ob%ships(n)%u /= missing_r ) then
162 apm_index = apm_index + 1
163 write(ounit,1000) ob%ships(n)%u, iv%ships(n)%u%inv, iv%ships(n)%u%error, &
164 'WNU', apm_plc, -888.88, iv%info(ships)%lat(1,n), iv%info(ships)%lon(1,n), &
165 (ob%ships(n)%p)/100., -888.88, apm_index
167 if ( iv%ships(n)%v%qc >= 0 .and. ob%ships(n)%v /= missing_r ) then
168 apm_index = apm_index + 1
169 write(ounit,1000) ob%ships(n)%v, iv%ships(n)%v%inv, iv%ships(n)%v%error, &
170 'WNV', apm_plc, -888.88, iv%info(ships)%lat(1,n), iv%info(ships)%lon(1,n), &
171 (ob%ships(n)%p)/100., -888.88, apm_index
173 if ( iv%ships(n)%t%qc >= 0 .and. ob%ships(n)%t /= missing_r ) then
174 apm_index = apm_index + 1
175 write(ounit,1000) ob%ships(n)%t, iv%ships(n)%t%inv, iv%ships(n)%t%error, &
176 'TMP', apm_plc, -888.88, iv%info(ships)%lat(1,n), iv%info(ships)%lon(1,n), &
177 (ob%ships(n)%p)/100., -888.88, apm_index
179 if ( iv%ships(n)%p%qc >= 0 .and. ob%ships(n)%p /= missing_r ) then
180 apm_index = apm_index + 1
181 write(ounit,1000) ob%ships(n)%p, iv%ships(n)%p%inv, iv%ships(n)%p%error, &
182 'PRS', apm_plc, -888.88, iv%info(ships)%lat(1,n), iv%info(ships)%lon(1,n), &
183 (ob%ships(n)%p)/100., -888.88, apm_index
185 if ( iv%ships(n)%q%qc >= 0 .and. ob%ships(n)%q /= missing_r ) then
186 apm_index = apm_index + 1
187 write(ounit,1000) ob%ships(n)%q, iv%ships(n)%q%inv, iv%ships(n)%q%error, &
188 'QVP',apm_plc, -888.88, iv%info(ships)%lat(1,n), iv%info(ships)%lon(1,n), &
189 (ob%ships(n)%p)/100., -888.88, apm_index
196 ! [4.1] Transfer Geo AMVs Obs:
198 if (iv%info(geoamv)%nlocal > 0) then
200 do n = 1, iv%info(geoamv)%nlocal
201 if (iv%info(geoamv)%proc_domain(1,n)) num_obs = num_obs + 1
203 if (num_obs > 0) then
205 do n = 1, iv%info(geoamv)%nlocal
206 if (iv%info(geoamv)%proc_domain(1,n)) then
207 num_obs = num_obs + 1
209 ! read(iv%info(geoamv)%platform(n)(4:6),*) apm_int
210 do k = 1, iv%info(geoamv)%levels(n)
211 if ( iv%geoamv(n)%u(k)%qc >= 0 .and. ob%geoamv(n)%u(k) /= missing_r ) then
212 apm_index = apm_index + 1
213 write(ounit,1000) ob%geoamv(n)%u(k), iv%geoamv(n)%u(k)%inv, iv%geoamv(n)%u(k)%error, &
214 'WNU', apm_plc, -888.88, iv%info(geoamv)%lat(1,n), iv%info(geoamv)%lon(1,n), &
215 (iv%geoamv(n)%p(k))/100., -888.88, apm_index
217 if ( iv%geoamv(n)%v(k)%qc >= 0 .and. ob%geoamv(n)%v(k) /= missing_r ) then
218 apm_index = apm_index + 1
219 write(ounit,1000) ob%geoamv(n)%v(k), iv%geoamv(n)%v(k)%inv, iv%geoamv(n)%v(k)%error, &
220 'WNV', apm_plc, -888.88, iv%info(geoamv)%lat(1,n), iv%info(geoamv)%lon(1,n), &
221 (iv%geoamv(n)%p(k))/100., -888.88, apm_index
229 ! [4.2] Transfer Polar AMVs Obs:
231 if (iv%info(polaramv)%nlocal > 0) then
233 do n = 1, iv%info(polaramv)%nlocal
234 if (iv%info(polaramv)%proc_domain(1,n)) num_obs = num_obs + 1
236 if (num_obs > 0) then
238 do n = 1, iv%info(polaramv)%nlocal
239 if (iv%info(polaramv)%proc_domain(1,n)) then
240 num_obs = num_obs + 1
242 ! read(iv%info(polaramv)%platform(n)(4:6),*) apm_int
243 do k = 1, iv%info(polaramv)%levels(n)
244 if ( iv%polaramv(n)%u(k)%qc >= 0 .and. ob%polaramv(n)%u(k) /= missing_r ) then
245 apm_index = apm_index + 1
246 write(ounit,1000) ob%polaramv(n)%u(k), iv%polaramv(n)%u(k)%inv, iv%polaramv(n)%u(k)%error, &
247 'WNU', apm_plc, -888.88, iv%info(polaramv)%lat(1,n), iv%info(polaramv)%lon(1,n), &
248 (iv%polaramv(n)%p(k))/100., -888.88, apm_index
250 if ( iv%polaramv(n)%v(k)%qc >= 0 .and. ob%polaramv(n)%v(k) /= missing_r ) then
251 apm_index = apm_index + 1
252 write(ounit,1000) ob%polaramv(n)%v(k), iv%polaramv(n)%v(k)%inv, iv%polaramv(n)%v(k)%error, &
253 'WNV', apm_plc, -888.88, iv%info(polaramv)%lat(1,n), iv%info(polaramv)%lon(1,n), &
254 (iv%polaramv(n)%p(k))/100., -888.88, apm_index
262 ! [5] Transfer gpspw obs:
264 if (iv%info(gpspw)%nlocal > 0) then
266 do n = 1, iv%info(gpspw)%nlocal
267 if (iv%info(gpspw)%proc_domain(1,n)) num_obs = num_obs + 1
269 if (num_obs > 0) then
271 do n = 1, iv%info(gpspw)%nlocal
272 if (iv%info(gpspw)%proc_domain(1,n)) then
273 num_obs = num_obs + 1
275 ! read(iv%info(gpspw)%platform(n)(4:6),*) apm_int
276 if ( iv%gpspw(n)%tpw%qc >= 0 .and. ob%gpspw(n)%tpw /= missing_r ) then
277 apm_index = apm_index + 1
278 write(ounit,1000) ob%gpspw(n)%tpw, iv%gpspw(n)%tpw%inv, iv%gpspw(n)%tpw%error, &
279 'PWT', apm_plc, -888.88, iv%info(gpspw)%lat(1,n), iv%info(gpspw)%lon(1,n), &
280 -888.88, -888.88, apm_index
287 ! [6] Transfer sonde obs:
289 if (iv%info(sound)%nlocal > 0) then
291 do n = 1, iv%info(sound)%nlocal
292 if (iv%info(sound)%proc_domain(1,n)) num_obs = num_obs + 1
294 if (num_obs > 0) then
296 do n = 1, iv%info(sound)%nlocal
297 if (iv%info(sound)%proc_domain(1,n)) then
298 num_obs = num_obs + 1
300 ! read(iv%info(sound)%platform(n)(4:6),*) apm_int
301 do k = 1, iv%info(sound)%levels(n)
302 if ( iv%sound(n)%u(k)%qc >= 0 .and. ob%sound(n)%u(k) /= missing_r ) then
303 apm_index = apm_index + 1
304 write(ounit,1000) ob%sound(n)%u(k), iv%sound(n)%u(k)%inv, iv%sound(n)%u(k)%error, &
305 'WNU', apm_plc, -888.88, iv%info(sound)%lat(1,n), iv%info(sound)%lon(1,n), &
306 (iv%sound(n)%p(k))/100., -888.88, apm_index
308 if ( iv%sound(n)%v(k)%qc >= 0 .and. ob%sound(n)%v(k) /= missing_r ) then
309 apm_index = apm_index + 1
310 write(ounit,1000) ob%sound(n)%v(k), iv%sound(n)%v(k)%inv, iv%sound(n)%v(k)%error, &
311 'WNV', apm_plc, -888.88, iv%info(sound)%lat(1,n), iv%info(sound)%lon(1,n), &
312 (iv%sound(n)%p(k))/100., -888.88, apm_index
314 if ( iv%sound(n)%t(k)%qc >= 0 .and. ob%sound(n)%t(k) /= missing_r ) then
315 apm_index = apm_index + 1
316 write(ounit,1000) ob%sound(n)%t(k), iv%sound(n)%t(k)%inv, iv%sound(n)%t(k)%error, &
317 'TMP', apm_plc, -888.88, iv%info(sound)%lat(1,n), iv%info(sound)%lon(1,n), &
318 (iv%sound(n)%p(k))/100., -888.88, apm_index
320 if ( iv%sound(n)%q(k)%qc >= 0 .and. ob%sound(n)%q(k) /= missing_r ) then
321 apm_index = apm_index + 1
322 write(ounit,1000) ob%sound(n)%q(k), iv%sound(n)%q(k)%inv, iv%sound(n)%q(k)%error, &
323 'QVP', apm_plc, -888.88, iv%info(sound)%lat(1,n), iv%info(sound)%lon(1,n), &
324 (iv%sound(n)%p(k))/100., -888.88, apm_index
332 if (iv%info(sonde_sfc)%nlocal > 0) then
334 do n = 1, iv%info(sonde_sfc)%nlocal
335 if (iv%info(sonde_sfc)%proc_domain(1,n)) num_obs = num_obs + 1
337 if (num_obs > 0) then
339 do n = 1, iv%info(sonde_sfc)%nlocal
340 if (iv%info(sonde_sfc)%proc_domain(1,n)) then
341 num_obs = num_obs + 1
343 ! read(iv%info(sonde_sfc)%platform(n)(4:6),*) apm_int
344 if ( iv%sonde_sfc(n)%u%qc >= 0 .and. ob%sonde_sfc(n)%u /= missing_r ) then
345 apm_index = apm_index + 1
346 write(ounit,1000) ob%sonde_sfc(n)%u, iv%sonde_sfc(n)%u%inv, iv%sonde_sfc(n)%u%error, &
347 'WNU', apm_plc, -888.88, iv%info(sonde_sfc)%lat(1,n), iv%info(sonde_sfc)%lon(1,n), &
348 (ob%sonde_sfc(n)%p)/100., -888.88, apm_index
350 if ( iv%sonde_sfc(n)%v%qc >= 0 .and. ob%sonde_sfc(n)%v /= missing_r ) then
351 apm_index = apm_index + 1
352 write(ounit,1000) ob%sonde_sfc(n)%v, iv%sonde_sfc(n)%v%inv, iv%sonde_sfc(n)%v%error, &
353 'WNV', apm_plc, -888.88, iv%info(sonde_sfc)%lat(1,n), iv%info(sonde_sfc)%lon(1,n), &
354 (ob%sonde_sfc(n)%p)/100., -888.88, apm_index
356 if ( iv%sonde_sfc(n)%t%qc >= 0 .and. ob%sonde_sfc(n)%t /= missing_r ) then
357 apm_index = apm_index + 1
358 write(ounit,1000) ob%sonde_sfc(n)%t, iv%sonde_sfc(n)%t%inv, iv%sonde_sfc(n)%t%error, &
359 'TMP', apm_plc, -888.88, iv%info(sonde_sfc)%lat(1,n), iv%info(sonde_sfc)%lon(1,n), &
360 (ob%sonde_sfc(n)%p)/100., -888.88, apm_index
362 if ( iv%sonde_sfc(n)%p%qc >= 0 .and. ob%sonde_sfc(n)%p /= missing_r ) then
363 apm_index = apm_index + 1
364 write(ounit,1000) ob%sonde_sfc(n)%p, iv%sonde_sfc(n)%p%inv, iv%sonde_sfc(n)%p%error, &
365 'PRS', apm_plc, -888.88, iv%info(sonde_sfc)%lat(1,n), iv%info(sonde_sfc)%lon(1,n), &
366 (ob%sonde_sfc(n)%p)/100., -888.88, apm_index
368 if ( iv%sonde_sfc(n)%q%qc >= 0 .and. ob%sonde_sfc(n)%q /= missing_r ) then
369 apm_index = apm_index + 1
370 write(ounit,1000) ob%sonde_sfc(n)%q, iv%sonde_sfc(n)%q%inv, iv%sonde_sfc(n)%q%error, &
371 'QVP', apm_plc, -888.88, iv%info(sonde_sfc)%lat(1,n), iv%info(sonde_sfc)%lon(1,n), &
372 (ob%sonde_sfc(n)%p)/100., -888.88, apm_index
379 ! [7] Transfer airep obs:
381 if (iv%info(airep)%nlocal > 0) then
383 do n = 1, iv%info(airep)%nlocal
384 if (iv%info(airep)%proc_domain(1,n)) num_obs = num_obs + 1
386 if (num_obs > 0) then
388 do n = 1, iv%info(airep)%nlocal
389 if (iv%info(airep)%proc_domain(1,n)) then
390 num_obs = num_obs + 1
392 ! read(iv%info(airep)%platform(n)(4:6),*) apm_int
393 do k = 1, iv%info(airep)%levels(n)
394 if ( iv%airep(n)%u(k)%qc >= 0 .and. ob%airep(n)%u(k) /= missing_r ) then
395 apm_index = apm_index + 1
396 write(ounit,1000) ob%airep(n)%u(k), iv%airep(n)%u(k)%inv, iv%airep(n)%u(k)%error, &
397 'WNU', apm_plc, -888.88, iv%info(airep)%lat(1,n), iv%info(airep)%lon(1,n), &
398 (iv%airep(n)%p(k))/100., -888.88, apm_index
400 if ( iv%airep(n)%v(k)%qc >= 0 .and. ob%airep(n)%v(k) /= missing_r ) then
401 apm_index = apm_index + 1
402 write(ounit,1000) ob%airep(n)%v(k), iv%airep(n)%v(k)%inv, iv%airep(n)%v(k)%error, &
403 'WNV', apm_plc, -888.88, iv%info(airep)%lat(1,n), iv%info(airep)%lon(1,n), &
404 (iv%airep(n)%p(k))/100., -888.88, apm_index
406 if ( iv%airep(n)%t(k)%qc >= 0 .and. ob%airep(n)%t(k) /= missing_r ) then
407 apm_index = apm_index + 1
408 write(ounit,1000) ob%airep(n)%t(k), iv%airep(n)%t(k)%inv, iv%airep(n)%t(k)%error, &
409 'TMP', apm_plc, -888.88, iv%info(airep)%lat(1,n), iv%info(airep)%lon(1,n), &
410 (iv%airep(n)%p(k))/100., -888.88, apm_index
418 ! [8] Transfer pilot obs:
420 if (iv%info(pilot)%nlocal > 0) then
422 do n = 1, iv%info(pilot)%nlocal
423 if (iv%info(pilot)%proc_domain(1,n)) num_obs = num_obs + 1
425 if (num_obs > 0) then
427 do n = 1, iv%info(pilot)%nlocal
428 if (iv%info(pilot)%proc_domain(1,n)) then
429 num_obs = num_obs + 1
431 ! read(iv%info(pilot)%platform(n)(4:6),*) apm_int
432 do k = 1, iv%info(pilot)%levels(n)
433 if ( iv%pilot(n)%u(k)%qc >= 0 .and. ob%pilot(n)%u(k) /= missing_r ) then
434 apm_index = apm_index + 1
435 write(ounit,1000) ob%pilot(n)%u(k), iv%pilot(n)%u(k)%inv, iv%pilot(n)%u(k)%error, &
436 'WNU', apm_plc, -888.88, iv%info(pilot)%lat(1,n), iv%info(pilot)%lon(1,n), &
437 (iv%pilot(n)%p(k))/100, -888.88, apm_index
439 if ( iv%pilot(n)%v(k)%qc >= 0 .and. ob%pilot(n)%v(k) /= missing_r ) then
440 apm_index = apm_index + 1
441 write(ounit,1000) ob%pilot(n)%v(k), iv%pilot(n)%v(k)%inv, iv%pilot(n)%v(k)%error, &
442 'WNV', apm_plc, -888.88, iv%info(pilot)%lat(1,n), iv%info(pilot)%lon(1,n), &
443 (iv%pilot(n)%p(k))/100., -888.88, apm_index
451 ! [9] Transfer SSM/I obs:SSMI:
453 if (iv%info(ssmi_rv)%nlocal > 0) then
455 do n = 1, iv%info(ssmi_rv)%nlocal
456 if (iv%info(ssmi_rv)%proc_domain(1,n)) num_obs = num_obs + 1
458 if (num_obs > 0) then
460 do n = 1, iv%info(ssmi_rv)%nlocal
461 if (iv%info(ssmi_rv)%proc_domain(1,n)) then
462 num_obs = num_obs + 1
464 ! read(iv%info(ssmi_rv)%platform(n)(4:6),*) apm_int
465 if ( iv%ssmi_rv(n)%speed%qc >= 0 .and. ob%ssmi_rv(n)%speed /= missing_r ) then
466 apm_index = apm_index + 1
467 write(ounit,1000) ob%ssmi_rv(n)%speed, iv%ssmi_rv(n)%speed%inv, &
468 iv%ssmi_rv(n)%speed%error, &
469 'SPD', apm_plc, -888.88, iv%info(ssmi_rv)%lat(1,n), iv%info(ssmi_rv)%lon(1,n), &
470 -888.88, -888.88, apm_index
472 if ( iv%ssmi_rv(n)%tpw%qc >= 0 .and. ob%ssmi_rv(n)%tpw /= missing_r ) then
473 apm_index = apm_index + 1
474 write(ounit,1000) ob%ssmi_rv(n)%tpw, iv%ssmi_rv(n)%tpw%inv, &
475 iv%ssmi_rv(n)%tpw%error, &
476 'PWT', apm_plc, -888.88, iv%info(ssmi_rv)%lat(1,n), iv%info(ssmi_rv)%lon(1,n), &
477 -888.88, -888.88, apm_index
484 ! SSM/I TB not coded.
486 ! [10] Transfer satem obs:
488 if (iv%info(satem)%nlocal > 0) then
490 do n = 1, iv%info(satem)%nlocal
491 if (iv%info(satem)%proc_domain(1,n)) num_obs = num_obs + 1
493 if (num_obs > 0) then
495 do n = 1, iv%info(satem)%nlocal
496 if (iv%info(satem)%proc_domain(1,n)) then
497 num_obs = num_obs + 1
499 ! read(iv%info(satem)%platform(n)(4:6),*) apm_int
500 do k = 1, iv%info(satem)%levels(n)
501 if ( iv%satem(n)%thickness(k)%qc >= 0 .and. ob%satem(n)%thickness(k) /= missing_r ) then
502 apm_index = apm_index + 1
503 write(ounit,1000) ob%satem(n)%thickness(k), iv%satem(n)%thickness(k)%inv, &
504 iv%satem(n)%thickness(k)%error, &
505 'TCK', apm_plc, -888.88, iv%info(satem)%lat(1,n), iv%info(satem)%lon(1,n), &
506 (iv%satem(n)%p(k))/100., -888.88, apm_index
514 ! SSMT1 SSMT2 not coded.
516 ! [11] Transfer scatterometer obs:
518 if (iv%info(qscat)%nlocal > 0) then
520 do n = 1, iv%info(qscat)%nlocal
521 if (iv%info(qscat)%proc_domain(1,n)) num_obs = num_obs + 1
523 if (num_obs > 0) then
525 do n = 1, iv%info(qscat)%nlocal
526 if (iv%info(qscat)%proc_domain(1,n)) then
527 num_obs = num_obs + 1
529 ! read(iv%info(qscat)%platform(n)(4:6),*) apm_int
530 if ( iv%qscat(n)%u%qc >= 0 .and. ob%qscat(n)%u /= missing_r ) then
531 apm_index = apm_index + 1
532 write(ounit,1000) ob%qscat(n)%u, iv%qscat(n)%u%inv, iv%qscat(n)%u%error, &
533 'WNU', apm_plc, -888.88, iv%info(qscat)%lat(1,n), iv%info(qscat)%lon(1,n), &
534 -888.88, -888.88, apm_index
536 if ( iv%qscat(n)%v%qc >= 0 .and. ob%qscat(n)%v /= missing_r ) then
537 apm_index = apm_index + 1
538 write(ounit,1000) ob%qscat(n)%v, iv%qscat(n)%v%inv, iv%qscat(n)%v%error, &
539 'WNV', apm_plc, -888.88, iv%info(qscat)%lat(1,n), iv%info(qscat)%lon(1,n), &
540 -888.88, -888.88, apm_index
547 ! [12] Transfer profiler obs:
549 if (iv%info(profiler)%nlocal > 0) then
551 do n = 1, iv%info(profiler)%nlocal
552 if (iv%info(profiler)%proc_domain(1,n)) num_obs = num_obs + 1
554 if (num_obs > 0) then
556 do n = 1, iv%info(profiler)%nlocal
557 if (iv%info(profiler)%proc_domain(1,n)) then
558 num_obs = num_obs + 1
560 ! read(iv%info(profiler)%platform(n)(4:6),*) apm_int
561 do k = 1, iv%info(profiler)%levels(n)
562 if ( iv%profiler(n)%u(k)%qc >= 0 .and. ob%profiler(n)%u(k) /= missing_r ) then
563 apm_index = apm_index + 1
564 write(ounit,1000) ob%profiler(n)%u(k), iv%profiler(n)%u(k)%inv, iv%profiler(n)%u(k)%error, &
565 'WNU', apm_plc, -888.88, iv%info(profiler)%lat(1,n), iv%info(profiler)%lon(1,n), &
566 (iv%profiler(n)%p(k))/100., -888.88, apm_index
568 if ( iv%profiler(n)%v(k)%qc >= 0 .and. ob%profiler(n)%v(k) /= missing_r ) then
569 apm_index = apm_index + 1
570 write(ounit,1000) ob%profiler(n)%v(k), iv%profiler(n)%v(k)%inv, iv%profiler(n)%v(k)%error, &
571 'WNV', apm_plc, -888.88, iv%info(profiler)%lat(1,n), iv%info(profiler)%lon(1,n), &
572 (iv%profiler(n)%p(k))/100., -888.88, apm_index
582 if (iv%info(buoy)%nlocal > 0) then
584 do n = 1, iv%info(buoy)%nlocal
585 if (iv%info(buoy)%proc_domain(1,n)) num_obs = num_obs + 1
587 if (num_obs > 0) then
589 do n = 1, iv%info(buoy)%nlocal
590 if (iv%info(buoy)%proc_domain(1,n)) then
591 num_obs = num_obs + 1
593 ! read(iv%info(buoy)%platform(n)(4:6),*) apm_int
594 if ( iv%buoy(n)%u%qc >= 0 .and. ob%buoy(n)%u /= missing_r ) then
595 apm_index = apm_index + 1
596 write(ounit,1000) ob%buoy(n)%u, iv%buoy(n)%u%inv, iv%buoy(n)%u%error, &
597 'WNU', apm_plc, -888.88, iv%info(buoy)%lat(1,n), iv%info(buoy)%lon(1,n), &
598 (ob%buoy(n)%p)/100., -888.88, apm_index
600 if ( iv%buoy(n)%v%qc >= 0 .and. ob%buoy(n)%v /= missing_r ) then
601 apm_index = apm_index + 1
602 write(ounit,1000) ob%buoy(n)%v, iv%buoy(n)%v%inv, iv%buoy(n)%v%error, &
603 'WNV', apm_plc, -888.88, iv%info(buoy)%lat(1,n), iv%info(buoy)%lon(1,n), &
604 (ob%buoy(n)%p)/100., -888.88, apm_index
606 if ( iv%buoy(n)%t%qc >= 0 .and. ob%buoy(n)%t /= missing_r ) then
607 apm_index = apm_index + 1
608 write(ounit,1000) ob%buoy(n)%t, iv%buoy(n)%t%inv, iv%buoy(n)%t%error, &
609 'TMP', apm_plc, -888.88, iv%info(buoy)%lat(1,n), iv%info(buoy)%lon(1,n), &
610 (ob%buoy(n)%p)/100., -888.88, apm_index
612 if ( iv%buoy(n)%p%qc >= 0 .and. ob%buoy(n)%p /= missing_r ) then
613 apm_index = apm_index + 1
614 write(ounit,1000) ob%buoy(n)%p, iv%buoy(n)%p%inv, iv%buoy(n)%p%error, &
615 'PRS', apm_plc, -888.88, iv%info(buoy)%lat(1,n), iv%info(buoy)%lon(1,n), &
616 (ob%buoy(n)%p)/100., -888.88, apm_index
618 if ( iv%buoy(n)%q%qc >= 0 .and. ob%buoy(n)%q /= missing_r ) then
619 apm_index = apm_index + 1
620 write(ounit,1000) ob%buoy(n)%q, iv%buoy(n)%q%inv, iv%buoy(n)%q%error, &
621 'QVP', apm_plc, -888.88, iv%info(buoy)%lat(1,n), iv%info(buoy)%lon(1,n), &
622 (ob%buoy(n)%p)/100., -888.88, apm_index
629 ! Transfer TC bogus obs:
631 if (iv%info(bogus)%nlocal > 0) then
633 do n = 1, iv%info(bogus)%nlocal
634 if (iv%info(bogus)%proc_domain(1,n)) num_obs = num_obs + 1
636 if (num_obs > 0) then
638 do n = 1, iv%info(bogus)%nlocal
639 if (iv%info(bogus)%proc_domain(1,n)) then
640 num_obs = num_obs + 1
642 ! read(iv%info(bogus)%platform(n)(4:6),*) apm_int
643 do k = 1, iv%info(bogus)%levels(n)
644 if ( iv%bogus(n)%u(k)%qc >= 0 .and. ob%bogus(n)%u(k) /= missing_r ) then
645 apm_index = apm_index + 1
646 write(ounit,1000) ob%bogus(n)%u(k), iv%bogus(n)%u(k)%inv, iv%bogus(n)%u(k)%error, &
647 'WNU', apm_plc, -888.88, iv%info(bogus)%lat(1,n), iv%info(bogus)%lon(1,n), &
648 (iv%bogus(n)%p(k))/100., -888.88, apm_index
650 if ( iv%bogus(n)%v(k)%qc >= 0 .and. ob%bogus(n)%v(k) /= missing_r ) then
651 apm_index = apm_index + 1
652 write(ounit,1000) ob%bogus(n)%v(k), iv%bogus(n)%v(k)%inv, iv%bogus(n)%v(k)%error, &
653 'WNV', apm_plc, -888.88, iv%info(bogus)%lat(1,n), iv%info(bogus)%lon(1,n), &
654 (iv%bogus(n)%p(k))/100., -888.88, apm_index
656 if ( iv%bogus(n)%t(k)%qc >= 0 .and. ob%bogus(n)%t(k) /= missing_r ) then
657 apm_index = apm_index + 1
658 write(ounit,1000) ob%bogus(n)%t(k), iv%bogus(n)%t(k)%inv, iv%bogus(n)%t(k)%error, &
659 'TMP', apm_plc, -888.88, iv%info(bogus)%lat(1,n), iv%info(bogus)%lon(1,n), &
660 (iv%bogus(n)%p(k))/100., -888.88, apm_index
662 if ( iv%bogus(n)%q(k)%qc >= 0 .and. ob%bogus(n)%q(k) /= missing_r ) then
663 apm_index = apm_index + 1
664 write(ounit,1000) ob%bogus(n)%q(k), iv%bogus(n)%q(k)%inv, iv%bogus(n)%q(k)%error, &
665 'QVP', apm_plc, -888.88, iv%info(bogus)%lat(1,n), iv%info(bogus)%lon(1,n), &
666 (iv%bogus(n)%p(k))/100., -888.88, apm_index
674 ! Transfer AIRS retrievals:
676 if (iv%info(airsr)%nlocal > 0) then
678 do n = 1, iv%info(airsr)%nlocal
679 if (iv%info(airsr)%proc_domain(1,n)) num_obs = num_obs + 1
681 if (num_obs > 0) then
683 do n = 1, iv%info(airsr)%nlocal
684 if (iv%info(airsr)%proc_domain(1,n)) then
685 num_obs = num_obs + 1
687 ! read(iv%info(airsr)%platform(n)(4:6),*) apm_int
688 do k = 1, iv%info(airsr)%levels(n)
689 if ( iv%airsr(n)%t(k)%qc >= 0 .and. ob%airsr(n)%t(k) /= missing_r ) then
690 apm_index = apm_index + 1
691 write(ounit,1000) ob%airsr(n)%t(k), iv%airsr(n)%t(k)%inv, iv%airsr(n)%t(k)%error, &
692 'TMP', apm_plc, -888.88, iv%info(airsr)%lat(1,n), iv%info(airsr)%lon(1,n), &
693 (iv%airsr(n)%p(k))/100., -888.88, apm_index
695 if ( iv%airsr(n)%q(k)%qc >= 0 .and. ob%airsr(n)%q(k) /= missing_r ) then
696 apm_index = apm_index + 1
697 write(ounit,1000) ob%airsr(n)%q(k), iv%airsr(n)%q(k)%inv, iv%airsr(n)%q(k)%error, &
698 'QVP', apm_plc, -888.88, iv%info(airsr)%lat(1,n), iv%info(airsr)%lon(1,n), &
699 (iv%airsr(n)%p(k))/100., -888.88, apm_index
707 ! Transfer gpsref obs:
709 if (iv%info(gpsref)%nlocal > 0) then
711 do n = 1, iv%info(gpsref)%nlocal
712 if (iv%info(gpsref)%proc_domain(1,n)) num_obs = num_obs + 1
714 if (num_obs > 0) then
716 do n = 1, iv%info(gpsref)%nlocal
717 if (iv%info(gpsref)%proc_domain(1,n)) then
718 num_obs = num_obs + 1
720 ! read(iv%info(gpsref)%platform(n)(4:6),*) apm_int
721 do k = 1, iv%info(gpsref)%levels(n)
722 if ( iv%gpsref(n)%ref(k)%qc >= 0 .and. ob%gpsref(n)%ref(k) /= missing_r ) then
723 apm_index = apm_index + 1
724 write(ounit,1000) ob%gpsref(n)%ref(k), iv%gpsref(n)%ref(k)%inv, iv%gpsref(n)%ref(k)%error, &
725 'REF', apm_plc, -888.88, iv%info(gpsref)%lat(1,n), iv%info(gpsref)%lon(1,n), &
726 -888.88, -888.88, apm_index
734 ! Transfer tamdar obs:
736 if (iv%info(tamdar)%nlocal > 0) then
738 do n = 1, iv%info(tamdar)%nlocal
739 if (iv%info(tamdar)%proc_domain(1,n)) num_obs = num_obs + 1
741 if (num_obs > 0) then
743 do n = 1, iv%info(tamdar)%nlocal
744 if (iv%info(tamdar)%proc_domain(1,n)) then
745 num_obs = num_obs + 1
747 ! read(iv%info(tamdar)%platform(n)(4:6),*) apm_int
748 do k = 1, iv%info(tamdar)%levels(n)
749 if ( iv%tamdar(n)%u(k)%qc >= 0 .and. ob%tamdar(n)%u(k) /= missing_r ) then
750 apm_index = apm_index + 1
751 write(ounit,1000) ob%tamdar(n)%u(k), iv%tamdar(n)%u(k)%inv, iv%tamdar(n)%u(k)%error, &
752 'WNU', apm_plc, -888.88, iv%info(tamdar)%lat(1,n), iv%info(tamdar)%lon(1,n), &
753 (iv%tamdar(n)%p(k))/100., -888.88, apm_index
755 if ( iv%tamdar(n)%v(k)%qc >= 0 .and. ob%tamdar(n)%v(k) /= missing_r ) then
756 apm_index = apm_index + 1
757 write(ounit,1000) ob%tamdar(n)%v(k), iv%tamdar(n)%v(k)%inv, iv%tamdar(n)%v(k)%error, &
758 'WNV', apm_plc, -888.88, iv%info(tamdar)%lat(1,n), iv%info(tamdar)%lon(1,n), &
759 (iv%tamdar(n)%p(k))/100., -888.88, apm_index
761 if ( iv%tamdar(n)%t(k)%qc >= 0 .and. ob%tamdar(n)%t(k) /= missing_r ) then
762 apm_index = apm_index + 1
763 write(ounit,1000) ob%tamdar(n)%t(k), iv%tamdar(n)%t(k)%inv, iv%tamdar(n)%t(k)%error, &
764 'TMP', apm_plc, -888.88, iv%info(tamdar)%lat(1,n), iv%info(tamdar)%lon(1,n), &
765 (iv%tamdar(n)%p(k))/100., -888.88, apm_index
767 if ( iv%tamdar(n)%q(k)%qc >= 0 .and. ob%tamdar(n)%q(k) /= missing_r ) then
768 apm_index = apm_index + 1
769 write(ounit,1000) ob%tamdar(n)%q(k), iv%tamdar(n)%q(k)%inv, iv%tamdar(n)%q(k)%error, &
770 'QVP', apm_plc, -888.88, iv%info(tamdar)%lat(1,n), iv%info(tamdar)%lon(1,n), &
771 (iv%tamdar(n)%p(k))/100., -888.88, apm_index
779 if (iv%info(tamdar_sfc)%nlocal > 0) then
781 do n = 1, iv%info(tamdar_sfc)%nlocal
782 if (iv%info(tamdar_sfc)%proc_domain(1,n)) num_obs = num_obs + 1
784 if (num_obs > 0) then
786 do n = 1, iv%info(tamdar_sfc)%nlocal
787 if (iv%info(tamdar_sfc)%proc_domain(1,n)) then
788 num_obs = num_obs + 1
790 ! read(iv%info(tamdar_sfc)%platform(n)(4:6),*) apm_int
791 if ( iv%tamdar_sfc(n)%u%qc >= 0 .and. ob%tamdar_sfc(n)%u /= missing_r ) then
792 apm_index = apm_index + 1
793 write(ounit,1000) ob%tamdar_sfc(n)%u, iv%tamdar_sfc(n)%u%inv, iv%tamdar_sfc(n)%u%error, &
794 'WNU', apm_plc, -888.88, iv%info(tamdar_sfc)%lat(1,n), iv%info(tamdar_sfc)%lon(1,n), &
795 (ob%tamdar_sfc(n)%p)/100., -888.88, apm_index
797 if ( iv%tamdar_sfc(n)%v%qc >= 0 .and. ob%tamdar_sfc(n)%v /= missing_r ) then
798 apm_index = apm_index + 1
799 write(ounit,1000) ob%tamdar_sfc(n)%v, iv%tamdar_sfc(n)%v%inv, iv%tamdar_sfc(n)%v%error, &
800 'WNV', apm_plc, -888.88, iv%info(tamdar_sfc)%lat(1,n), iv%info(tamdar_sfc)%lon(1,n), &
801 (ob%tamdar_sfc(n)%p)/100., -888.88, apm_index
803 if ( iv%tamdar_sfc(n)%t%qc >= 0 .and. ob%tamdar_sfc(n)%t /= missing_r ) then
804 apm_index = apm_index + 1
805 write(ounit,1000) ob%tamdar_sfc(n)%t, iv%tamdar_sfc(n)%t%inv, iv%tamdar_sfc(n)%t%error, &
806 'TMP', apm_plc, -888.88, iv%info(tamdar_sfc)%lat(1,n), iv%info(tamdar_sfc)%lon(1,n), &
807 (ob%tamdar_sfc(n)%p)/100., -888.88, apm_index
809 if ( iv%tamdar_sfc(n)%p%qc >= 0 .and. ob%tamdar_sfc(n)%p /= missing_r ) then
810 apm_index = apm_index + 1
811 write(ounit,1000) ob%tamdar_sfc(n)%p, iv%tamdar_sfc(n)%p%inv, iv%tamdar_sfc(n)%p%error, &
812 'PRS', apm_plc, -888.88, iv%info(tamdar_sfc)%lat(1,n), iv%info(tamdar_sfc)%lon(1,n), &
813 (ob%tamdar_sfc(n)%p)/100., -888.88, apm_index
815 if ( iv%tamdar_sfc(n)%q%qc >= 0 .and. ob%tamdar_sfc(n)%q /= missing_r ) then
816 apm_index = apm_index + 1
817 write(ounit,1000) ob%tamdar_sfc(n)%q, iv%tamdar_sfc(n)%q%inv, iv%tamdar_sfc(n)%q%error, &
818 'QVP', apm_plc, -888.88, iv%info(tamdar_sfc)%lat(1,n), iv%info(tamdar_sfc)%lon(1,n), &
819 (ob%tamdar_sfc(n)%p)/100., -888.88, apm_index
827 call da_free_unit(ounit)
829 if (trace_use) call da_trace_exit("da_write_obs_etkf")
831 end subroutine da_write_obs_etkf