wrf svn trunk commit r4103
[wrffire.git] / wrfv2_fire / var / da / da_obs / da_fill_obs_structures.inc
blob9c74f272aa845b114cc7df10b9c158b9d773d518
1 subroutine da_fill_obs_structures(iv, ob, uvq_direct)
3    !----------------------------------------------------------------------------   
4    ! Purpose: Allocates observation structure and fills it from iv.
5    !----------------------------------------------------------------------------   
7    implicit none
9    type (iv_type), intent(inout) :: iv   ! Obs and header structure.
10    type (y_type), intent(out)    :: ob   ! (Smaller) observation structure.
11    logical, optional             :: uvq_direct  !flag for having direct u,v,q obs
13    integer :: n, k     ! Loop counters.
14    real    :: rh_error ! RH obs. error.
15    real    :: q_error  ! q obs. error.
16    integer :: i,j
17    logical :: outside
19    if (trace_use) call da_trace_entry("da_fill_obs_structures")
21    !---------------------------------------------------------------------------
22    ! Initialise obs error factors (which will be overwritten in use_obs_errfac)
23    !---------------------------------------------------------------------------
25    iv % synop_ef_u = 1.0
26    iv % synop_ef_v = 1.0
27    iv % synop_ef_t = 1.0
28    iv % synop_ef_p = 1.0
29    iv % synop_ef_q = 1.0
31    iv % metar_ef_u = 1.0
32    iv % metar_ef_v = 1.0
33    iv % metar_ef_t = 1.0
34    iv % metar_ef_p = 1.0
35    iv % metar_ef_q = 1.0
37    iv % ships_ef_u = 1.0
38    iv % ships_ef_v = 1.0
39    iv % ships_ef_t = 1.0
40    iv % ships_ef_p = 1.0
41    iv % ships_ef_q = 1.0
43    iv % geoamv_ef_u = 1.0
44    iv % geoamv_ef_v = 1.0
46    iv % polaramv_ef_u = 1.0
47    iv % polaramv_ef_v = 1.0
49    iv % gpspw_ef_tpw = 1.0
51    iv % gpsref_ef_ref = 1.0
52    iv % gpsref_ef_p = 1.0
53    iv % gpsref_ef_t = 1.0
54    iv % gpsref_ef_q = 1.0
56    iv % sound_ef_u = 1.0
57    iv % sound_ef_v = 1.0
58    iv % sound_ef_t = 1.0
59    iv % sound_ef_q = 1.0
61    iv % airep_ef_u = 1.0
62    iv % airep_ef_v = 1.0
63    iv % airep_ef_t = 1.0
65    iv % pilot_ef_u = 1.0
66    iv % pilot_ef_v = 1.0
68    iv % ssmir_ef_speed = 1.0
69    iv % ssmir_ef_tpw = 1.0
71    iv % satem_ef_thickness = 1.0
73    iv % ssmt1_ef_t = 1.0
75    iv % ssmt2_ef_rh = 1.0
77    iv % qscat_ef_u = 1.0
78    iv % qscat_ef_v = 1.0
80    iv % profiler_ef_u = 1.0
81    iv % profiler_ef_v = 1.0
82    
83    iv % buoy_ef_u = 1.0
84    iv % buoy_ef_v = 1.0
85    iv % buoy_ef_t = 1.0
86    iv % buoy_ef_p = 1.0
87    iv % buoy_ef_q = 1.0
89    iv % radar_ef_rv = 1.0
90    iv % radar_ef_rf = 1.0
92    iv % bogus_ef_u = 1.0
93    iv % bogus_ef_v = 1.0
94    iv % bogus_ef_t = 1.0
95    iv % bogus_ef_p = 1.0
96    iv % bogus_ef_q = 1.0
97    iv % bogus_ef_slp = 1.0
99    iv % airsr_ef_t = 1.0
100    iv % airsr_ef_q = 1.0
102    iv % mtgirs_ef_u = 1.0
103    iv % mtgirs_ef_v = 1.0
104    iv % mtgirs_ef_t = 1.0
105    iv % mtgirs_ef_q = 1.0
107    iv % tamdar_ef_u = 1.0
108    iv % tamdar_ef_v = 1.0
109    iv % tamdar_ef_t = 1.0
110    iv % tamdar_ef_q = 1.0
112    iv % tamdar_sfc_ef_u = 1.0
113    iv % tamdar_sfc_ef_v = 1.0
114    iv % tamdar_sfc_ef_t = 1.0
115    iv % tamdar_sfc_ef_p = 1.0
116    iv % tamdar_sfc_ef_q = 1.0
118    !----------------------------------------------------------------------
119    ! [1.0] Allocate innovation vector and observation structures:
120    !----------------------------------------------------------------------
121    call da_allocate_y(iv, ob)
123    !----------------------------------------------------------------------
124    ! [2.0] Transfer observations:
125    !----------------------------------------------------------------------
127    ! [2.1] Transfer surface obs:
129    if (iv%info(synop)%nlocal > 0) then
130       do n = 1, iv%info(synop)%nlocal
131          ob % synop(n) % u = iv % synop(n) % u % inv
132          ob % synop(n) % v = iv % synop(n) % v % inv
133          ob % synop(n) % t = iv % synop(n) % t % inv
134          ob % synop(n) % q = iv % synop(n) % q % inv
135          ob % synop(n) % p = iv % synop(n) % p % inv
137          ! Calculate q error from rh error:
139          if (.not. present(uvq_direct) .or. (present(uvq_direct) .and. (.not. uvq_direct))) then   !cys_add
140          rh_error = iv%synop(n)%q%error ! q error is rh at this stage!
142          ! if((ob % synop(n) % p > iv%ptop) .AND. &
143          !    (ob % synop(n) % t > 100.0) .AND. &
144          !    (ob % synop(n) % q > 0.0) .AND. &
145          !    (iv % synop(n) % p % qc >= obs_qc_pointer) .and. &
146          !    (iv % synop(n) % t % qc >= obs_qc_pointer) .and. &
147          !    (iv % synop(n) % q % qc >= obs_qc_pointer)) then
148          call da_get_q_error(ob % synop(n) % p, &
149                               ob % synop(n) % t, &
150                               ob % synop(n) % q, &
151                               iv % synop(n) % t % error, &
152                               rh_error, iv % synop(n) % q % error)
153          if (iv%synop(n)% q % error == missing_r) iv%synop(n)% q % qc = missing_data
155          ! end if
156          end if   !cys_add
157       end do      
158    end if
160    ! [2.2] Transfer metar obs:
162    if (iv%info(metar)%nlocal > 0) then
163       do n = 1, iv%info(metar)%nlocal
164          ob % metar(n) % u = iv % metar(n) % u % inv
165          ob % metar(n) % v = iv % metar(n) % v % inv
166          ob % metar(n) % t = iv % metar(n) % t % inv
167          ob % metar(n) % q = iv % metar(n) % q % inv
168          ob % metar(n) % p = iv % metar(n) % p % inv
170          ! Calculate q error from rh error:
172          if (.not. present(uvq_direct) .or. (present(uvq_direct) .and. (.not. uvq_direct))) then   !cys_add
173          rh_error = iv%metar(n)%q%error ! q error is rh at this stage!
174          call da_get_q_error(iv % metar(n) % p % inv, &
175                               ob % metar(n) % t, &
176                               ob % metar(n) % q, &
177                               iv % metar(n) % t % error, &
178                               rh_error, q_error)
179          iv % metar(n) % q % error = q_error
180          if (iv%metar(n)% q % error == missing_r) &
181             iv%metar(n)% q % qc = missing_data
182          end if   !cys_add
183       end do
184    end if
186    ! [2.2] Transfer ships obs:
188    if (iv%info(ships)%nlocal > 0) then   
189       do n = 1, iv%info(ships)%nlocal
190          ob % ships(n) % u = iv % ships(n) % u % inv
191          ob % ships(n) % v = iv % ships(n) % v % inv
192          ob % ships(n) % t = iv % ships(n) % t % inv
193          ob % ships(n) % q = iv % ships(n) % q % inv
194          ob % ships(n) % p = iv % ships(n) % p % inv
196          ! Calculate q error from rh error:
198          if (.not. present(uvq_direct) .or. (present(uvq_direct) .and. (.not. uvq_direct))) then   !cys_add
199          rh_error = iv%ships(n)%q%error ! q error is rh at this stage!
200          call da_get_q_error(iv % ships(n) % p % inv, &
201                               ob % ships(n) % t, &
202                               ob % ships(n) % q, &
203                               iv % ships(n) % t % error, &
204                               rh_error, q_error)
205          iv % ships(n) % q % error = q_error
207          if(iv%ships(n)% q % error == missing_r) iv%ships(n)% q % qc = missing_data
208          end if   !cys_add
209       end do
210       
211    end if
213    ! [2.4.1] Transfer Geo. AMVs Obs:
215    if (iv%info(geoamv)%nlocal > 0) then
216       do n = 1, iv%info(geoamv)%nlocal
217          do k = 1, iv%info(geoamv)%levels(n)
218             ob % geoamv(n) % u(k) = iv % geoamv(n) % u(k) % inv
219             ob % geoamv(n) % v(k) = iv % geoamv(n) % v(k) % inv
220          end do
221       end do
222    end if
224    ! [2.4.2] Transfer  Polar AMVs Obs:
226    if (iv%info(polaramv)%nlocal > 0) then
227       do n = 1, iv%info(polaramv)%nlocal
228          do k = 1,iv%info(polaramv)%levels(n)
229             ob % polaramv(n) % u(k) = iv % polaramv(n) % u(k) % inv
230             ob % polaramv(n) % v(k) = iv % polaramv(n) % v(k) % inv
231          end do
232       end do
233    end if
235    ! [2.5] Transfer gpspw obs:
237    if (iv%info(gpspw)%nlocal > 0) then
238       do n = 1, iv%info(gpspw)%nlocal
239          ob % gpspw(n) % tpw = iv % gpspw(n) % tpw % inv
240       end do
242    end if
244    ! [2.6] Transfer GPS REF obs:
246    if (iv%info(gpsref)%nlocal > 0) then
247       do n = 1, iv%info(gpsref)%nlocal
248          do k = 1, iv%info(gpsref)%levels(n)
249             ob % gpsref(n) % ref(k) = iv % gpsref(n) % ref(k) % inv
250             ob % gpsref(n) %   p(k) = iv % gpsref(n) %   p(k) % inv
251             ob % gpsref(n) %   t(k) = iv % gpsref(n) %   t(k) % inv
252             ob % gpsref(n) %   q(k) = iv % gpsref(n) %   q(k) % inv
253          end do
254       end do
255    end if
257    ! [2.7] Transfer sonde obs:
259    if (iv%info(sound)%nlocal > 0) then
260       do n = 1, iv%info(sound)%nlocal
261          do k = 1, iv%info(sound)%levels(n)
262             ob % sound(n) % u(k) = iv % sound(n) % u(k) % inv
263             ob % sound(n) % v(k) = iv % sound(n) % v(k) % inv
264             ob % sound(n) % t(k) = iv % sound(n) % t(k) % inv
265             ob % sound(n) % q(k) = iv % sound(n) % q(k) % inv
267             ! Calculate q error from rh error:
269          if (.not. present(uvq_direct) .or. (present(uvq_direct) .and. (.not. uvq_direct))) then   !cys_add
270             rh_error = iv%sound(n)%q(k)%error ! q error is rh at this stage!
271             call da_get_q_error(iv % sound(n) % p(k), &
272                                  ob % sound(n) % t(k), &
273                                  ob % sound(n) % q(k), &
274                                  iv % sound(n) % t(k) % error, &
275                                  rh_error, q_error)
277             iv % sound(n) % q(k) % error = q_error
278          if (iv%sound(n)% q(k) % error == missing_r) &
279             iv%sound(n)% q(k) % qc = missing_data
280          end if   !cys_add
281          end do
282          ob % sonde_sfc(n) % u = iv % sonde_sfc(n) % u % inv
283          ob % sonde_sfc(n) % v = iv % sonde_sfc(n) % v % inv
284          ob % sonde_sfc(n) % t = iv % sonde_sfc(n) % t % inv
285          ob % sonde_sfc(n) % q = iv % sonde_sfc(n) % q % inv
286          ob % sonde_sfc(n) % p = iv % sonde_sfc(n) % p % inv
288          ! Calculate q error from rh error:
290          if (.not. present(uvq_direct) .or. (present(uvq_direct) .and. (.not. uvq_direct))) then   !cys_add
291          rh_error = iv%sonde_sfc(n)%q%error ! q error is rh at this stage!
292          call da_get_q_error(iv % sonde_sfc(n) % p % inv, &
293                               ob % sonde_sfc(n) % t, &
294                               ob % sonde_sfc(n) % q, &
295                               iv % sonde_sfc(n) % t % error, &
296                               rh_error, iv % sonde_sfc(n) % q % error)
297          if (iv%sonde_sfc(n)% q % error == missing_r) &
298             iv%sonde_sfc(n)% q % qc = missing_data
299          end if   !cys_add
300       end do
301    end if
303    ! [2.8] Transfer airep obs:
305    if (iv%info(airep)%nlocal > 0) then
306       do n = 1, iv%info(airep)%nlocal
307          do k = 1, iv%info(airep)%levels(n)
308             ob % airep(n) % u(k) = iv % airep(n) % u(k) % inv
309             ob % airep(n) % v(k) = iv % airep(n) % v(k) % inv
310             ob % airep(n) % t(k) = iv % airep(n) % t(k) % inv
311          end do
312       end do
313    end if
315    ! [2.9] Transfer pilot obs:
317    if (iv%info(pilot)%nlocal > 0) then
318       do n = 1, iv%info(pilot)%nlocal
319          do k = 1, iv%info(pilot)%levels(n)
320             ob % pilot(n) % u(k) = iv % pilot(n) % u(k) % inv
321             ob % pilot(n) % v(k) = iv % pilot(n) % v(k) % inv
322          end do
323       end do
324    end if
326    ! [2.10] Transfer SSM/I obs:SSMI:
328    if (iv%info(ssmi_rv)%nlocal > 0) then
329       do n = 1, iv%info(ssmi_rv)%nlocal
330          ob % ssmi_rv(n) % speed = iv % ssmi_rv(n) % speed % inv
331          ob % ssmi_rv(n) % tpw   = iv % ssmi_rv(n) % tpw % inv
332       end do
333    end if
335    if (iv%info(ssmi_tb)%nlocal > 0) then
336       do n = 1, iv%info(ssmi_tb)%nlocal
337          ob % ssmi_tb(n) % tb19v = iv % ssmi_tb(n) % tb19v % inv
338          ob % ssmi_tb(n) % tb19h = iv % ssmi_tb(n) % tb19h % inv
339          ob % ssmi_tb(n) % tb22v = iv % ssmi_tb(n) % tb22v % inv
340          ob % ssmi_tb(n) % tb37v = iv % ssmi_tb(n) % tb37v % inv
341          ob % ssmi_tb(n) % tb37h = iv % ssmi_tb(n) % tb37h % inv
342          ob % ssmi_tb(n) % tb85v = iv % ssmi_tb(n) % tb85v % inv
343          ob % ssmi_tb(n) % tb85h = iv % ssmi_tb(n) % tb85h % inv
344       end do
345    end if
347    ! [2.11] Transfer satem obs:
349    if (iv%info(satem)%nlocal > 0) then
350       do n = 1, iv%info(satem)%nlocal
351          do k = 1, iv%info(satem)%levels(n)
352             ob % satem(n) % thickness(k) = iv % satem(n) % thickness(k) % inv
353          end do
354       end do
355    end if
356    
357    ! [2.12] Transfer ssmt1 obs:
359    if (iv%info(ssmt1)%nlocal > 0) then
360       do n = 1, iv%info(ssmt1)%nlocal
361          do k = 1, iv%info(ssmt1)%levels(n)
362             ob % ssmt1(n) % t(k) = iv % ssmt1(n) % t(k) % inv
363          end do
364       end do
366    end if
368    ! [2.13] Transfer ssmt2 obs:
370    if (iv%info(ssmt2)%nlocal > 0) then
371       do n = 1, iv%info(ssmt2)%nlocal
372          do k = 1, iv%info(ssmt2)%levels(n)
373             ob % ssmt2(n) % rh(k) = iv % ssmt2(n) % rh(k) % inv
374          end do
375       end do
376    end if
377    
378    ! [2.14] Setup pseudo observations:
380    if (iv%info(pseudo)%nlocal > 0) call da_setup_pseudo_obs(iv, ob)
382    ! [2.15] Transfer scatterometer obs:
384    if (iv%info(qscat)%nlocal > 0) then
385       do n = 1, iv%info(qscat)%nlocal
386          ob % qscat(n) % u = iv % qscat(n) % u % inv
387          ob % qscat(n) % v = iv % qscat(n) % v % inv
388       end do     
389    end if
391    ! [2.16] Transfer profiler obs:
393    if (iv%info(profiler)%nlocal > 0) then
394       do n = 1, iv%info(profiler)%nlocal
395          do k = 1, iv%info(profiler)%levels(n)
396             ob % profiler(n) % u(k) = iv % profiler(n) % u(k) % inv
397             ob % profiler(n) % v(k) = iv % profiler(n) % v(k) % inv
398          end do
399       end do
400    end if
402    ! [2.17] Transfer buoy obs:
404    if (iv%info(buoy)%nlocal > 0) then
405       do n = 1, iv%info(buoy)%nlocal
406          ob % buoy(n) % p = iv % buoy(n) % p % inv
407       end do
408       do n = 1, iv%info(buoy)%nlocal
409          ob % buoy(n) % u = iv % buoy(n) % u % inv
410          ob % buoy(n) % v = iv % buoy(n) % v % inv
411          ob % buoy(n) % t = iv % buoy(n) % t % inv
412          ob % buoy(n) % q = iv % buoy(n) % q % inv
414          ! Calculate q error from rh error:
416          if (.not. present(uvq_direct) .or. (present(uvq_direct) .and. (.not. uvq_direct))) then   !cys_add
417          rh_error = iv%buoy(n)%q%error ! q error is rh at this stage!
418          call da_get_q_error(iv % buoy(n) % p % inv, &
419                               ob % buoy(n) % t, &
420                               ob % buoy(n) % q, &
421                               iv % buoy(n) % t % error, &
422                               rh_error, q_error)
423          iv % buoy(n) % q % error = q_error
425          if(iv%buoy (n)% q % error == missing_r) iv%buoy (n)% q % qc = missing_data
426          end if   !cys_add
427       end do
428    end if
430    ! [2.18] Transfer radar obs:
432    if (iv%info(radar)%nlocal > 0) then
433       do n = 1, iv%info(radar)%nlocal
434          do k = 1, iv%info(radar)%levels(n)
435             ! Copy observation variables:
436             ob % radar(n) % rv(k) = iv % radar(n) % rv(k) % inv
437            ob % radar(n) % rf(k) = iv % radar(n) % rf(k) % inv
438          end do
439       end do
440    end if
442    ! [2.19] Transfer TC bogus:
444    if (iv%info(bogus)%nlocal > 0) then
445       do n = 1, iv%info(bogus)%nlocal
446          do k = 1, iv%info(bogus)%levels(n)
448             ! Copy observation variables:
450             ob % bogus(n) % u(k) = iv % bogus(n) % u(k) % inv
451             ob % bogus(n) % v(k) = iv % bogus(n) % v(k) % inv
452             ob % bogus(n) % t(k) = iv % bogus(n) % t(k) % inv
453             ob % bogus(n) % q(k) = iv % bogus(n) % q(k) % inv
455             ! Calculate q error from rh error:
457             rh_error = iv%bogus(n)%q(k)%error ! q error is rh at this stage!
458             call da_get_q_error(iv % bogus(n) % p(k), &
459                                  ob % bogus(n) % t(k), &
460                                  ob % bogus(n) % q(k), &
461                                  iv % bogus(n) % t(k) % error, &
462                                  rh_error, q_error)
464             iv % bogus(n) % q(k) % error = q_error
465             if (iv%bogus(n)% q(k) % error == missing_r) &
466                iv%bogus(n)% q(k) % qc = missing_data
467          end do
468          ob % bogus(n) % slp = iv % bogus(n) % slp % inv
469       end do
470    end if
472    ! Transfer AIRS  retrievals:
474    if (iv%info(airsr)%nlocal > 0) then
475       do n = 1, iv%info(airsr)%nlocal
476          do k = 1, iv%info(airsr)%levels(n)
478             ! Copy observation variables:
480             ob % airsr(n) % t(k) = iv % airsr(n) % t(k) % inv
481             ob % airsr(n) % q(k) = iv % airsr(n) % q(k) % inv
483             ! Calculate q error from rh error:
485          if (.not. present(uvq_direct) .or. (present(uvq_direct) .and. (.not. uvq_direct))) then   !cys_add
486             rh_error = iv%airsr(n)%q(k)%error ! q error is rh at this stage!
487             call da_get_q_error(iv % airsr(n) % p(k), &
488                                  ob % airsr(n) % t(k), &
489                                  ob % airsr(n) % q(k), &
490                                  iv % airsr(n) % t(k) % error, &
491                                  rh_error, q_error)
493             iv % airsr(n) % q(k) % error = q_error
494             if (iv%airsr(n)% q(k) % error == missing_r) &
495                iv%airsr(n)% q(k) % qc = missing_data
496          end if   !cys_add
497          end do
498       end do
499    end if
500    if (iv%info(mtgirs)%nlocal > 0) then
501       do n = 1, iv%info(mtgirs)%nlocal
502          do k = 1, iv%info(mtgirs)%levels(n)
503             ob % mtgirs(n) % u(k) = iv % mtgirs(n) % u(k) % inv
504             ob % mtgirs(n) % v(k) = iv % mtgirs(n) % v(k) % inv
505             ob % mtgirs(n) % t(k) = iv % mtgirs(n) % t(k) % inv
506             ob % mtgirs(n) % q(k) = iv % mtgirs(n) % q(k) % inv
507          if (iv%mtgirs(n)% q(k) % error == missing_r) &
508             iv%mtgirs(n)% q(k) % qc = missing_data
509          end do
510       end do
511    end if
513    if (iv%info(tamdar)%nlocal > 0) then
514       do n = 1, iv%info(tamdar)%nlocal
515          do k = 1, iv%info(tamdar)%levels(n)
516             ob % tamdar(n) % u(k) = iv % tamdar(n) % u(k) % inv
517             ob % tamdar(n) % v(k) = iv % tamdar(n) % v(k) % inv
518             ob % tamdar(n) % t(k) = iv % tamdar(n) % t(k) % inv
519             ob % tamdar(n) % q(k) = iv % tamdar(n) % q(k) % inv
521          if (iv%tamdar(n)% u(k) % error == missing_r) &
522             iv%tamdar(n)% u(k) % qc = missing_data
523          if (iv%tamdar(n)% v(k) % error == missing_r) &
524             iv%tamdar(n)% v(k) % qc = missing_data
525          if (iv%tamdar(n)% t(k) % error == missing_r) &
526             iv%tamdar(n)% t(k) % qc = missing_data
528             ! Calculate q error from rh error:
530             rh_error = iv%tamdar(n)%q(k)%error ! q error is rh at this stage!
531             call da_get_q_error(iv % tamdar(n) % p(k), &
532                                 ob % tamdar(n) % t(k), &
533                                 ob % tamdar(n) % q(k), &
534                                 iv % tamdar(n) % t(k) % error, &
535                                 rh_error, q_error)
537             iv % tamdar(n) % q(k) % error = q_error
539          if (iv%tamdar(n)% q(k) % error == missing_r) &
540             iv%tamdar(n)% q(k) % qc = missing_data
541          end do
542          ob % tamdar_sfc(n) % u = iv % tamdar_sfc(n) % u % inv
543          ob % tamdar_sfc(n) % v = iv % tamdar_sfc(n) % v % inv
544          ob % tamdar_sfc(n) % t = iv % tamdar_sfc(n) % t % inv
545          ob % tamdar_sfc(n) % q = iv % tamdar_sfc(n) % q % inv
546          ob % tamdar_sfc(n) % p = iv % tamdar_sfc(n) % p % inv
548          if (iv%tamdar_sfc(n)% u % error == missing_r) &
549             iv%tamdar_sfc(n)% u % qc = missing_data
550          if (iv%tamdar_sfc(n)% v % error == missing_r) &
551             iv%tamdar_sfc(n)% v % qc = missing_data
552          if (iv%tamdar_sfc(n)% t % error == missing_r) &
553             iv%tamdar_sfc(n)% t % qc = missing_data
555          ! Calculate q error from rh error:
557          rh_error = iv%tamdar_sfc(n)%q%error ! q error is rh at this stage!
558          call da_get_q_error(iv % tamdar_sfc(n) % p % inv, &
559                               ob % tamdar_sfc(n) % t, &
560                               ob % tamdar_sfc(n) % q, &
561                               iv % tamdar_sfc(n) % t % error, &
562                               rh_error, iv % tamdar_sfc(n) % q % error)
563          if (iv%tamdar_sfc(n)% q % error == missing_r) &
564             iv%tamdar_sfc(n)% q % qc = missing_data
565       end do
566    end if
568    if (trace_use) call da_trace_exit("da_fill_obs_structures")
570 end subroutine da_fill_obs_structures