Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_obs / da_fill_obs_structures.inc
blob48a877b889a14ecb39ad523999f7810ee13e7952
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    real    :: geometric_h, geopotential_h
17    integer :: i,j
18    logical :: outside
20    if (trace_use) call da_trace_entry("da_fill_obs_structures")
22    !---------------------------------------------------------------------------
23    ! Initialise obs error factors (which will be overwritten in use_obs_errfac)
24    !---------------------------------------------------------------------------
26    iv % synop_ef_u = 1.0
27    iv % synop_ef_v = 1.0
28    iv % synop_ef_t = 1.0
29    iv % synop_ef_p = 1.0
30    iv % synop_ef_q = 1.0
32    iv % metar_ef_u = 1.0
33    iv % metar_ef_v = 1.0
34    iv % metar_ef_t = 1.0
35    iv % metar_ef_p = 1.0
36    iv % metar_ef_q = 1.0
38    iv % ships_ef_u = 1.0
39    iv % ships_ef_v = 1.0
40    iv % ships_ef_t = 1.0
41    iv % ships_ef_p = 1.0
42    iv % ships_ef_q = 1.0
44    iv % geoamv_ef_u = 1.0
45    iv % geoamv_ef_v = 1.0
47    iv % polaramv_ef_u = 1.0
48    iv % polaramv_ef_v = 1.0
50    iv % gpspw_ef_tpw = 1.0
52    iv % gpsref_ef_ref = 1.0
53    iv % gpsref_ef_p = 1.0
54    iv % gpsref_ef_t = 1.0
55    iv % gpsref_ef_q = 1.0
57    iv % gpseph_ef_eph = 1.0
59    iv % sound_ef_u = 1.0
60    iv % sound_ef_v = 1.0
61    iv % sound_ef_t = 1.0
62    iv % sound_ef_q = 1.0
64    iv % airep_ef_u = 1.0
65    iv % airep_ef_v = 1.0
66    iv % airep_ef_t = 1.0
67    iv % airep_ef_q = 1.0
69    iv % pilot_ef_u = 1.0
70    iv % pilot_ef_v = 1.0
72    iv % ssmir_ef_speed = 1.0
73    iv % ssmir_ef_tpw = 1.0
75    iv % satem_ef_thickness = 1.0
77    iv % ssmt1_ef_t = 1.0
79    iv % ssmt2_ef_rh = 1.0
81    iv % qscat_ef_u = 1.0
82    iv % qscat_ef_v = 1.0
84    iv % profiler_ef_u = 1.0
85    iv % profiler_ef_v = 1.0
86    
87    iv % buoy_ef_u = 1.0
88    iv % buoy_ef_v = 1.0
89    iv % buoy_ef_t = 1.0
90    iv % buoy_ef_p = 1.0
91    iv % buoy_ef_q = 1.0
93    iv % radar_ef_rv = 1.0
94    iv % radar_ef_rf = 1.0
96    iv % rain_ef_r  = 1.0
98    iv % bogus_ef_u = 1.0
99    iv % bogus_ef_v = 1.0
100    iv % bogus_ef_t = 1.0
101    iv % bogus_ef_p = 1.0
102    iv % bogus_ef_q = 1.0
103    iv % bogus_ef_slp = 1.0
105    iv % airsr_ef_t = 1.0
106    iv % airsr_ef_q = 1.0
108    iv % mtgirs_ef_u = 1.0
109    iv % mtgirs_ef_v = 1.0
110    iv % mtgirs_ef_t = 1.0
111    iv % mtgirs_ef_q = 1.0
113    iv % tamdar_ef_u = 1.0
114    iv % tamdar_ef_v = 1.0
115    iv % tamdar_ef_t = 1.0
116    iv % tamdar_ef_q = 1.0
118    iv % tamdar_sfc_ef_u = 1.0
119    iv % tamdar_sfc_ef_v = 1.0
120    iv % tamdar_sfc_ef_t = 1.0
121    iv % tamdar_sfc_ef_p = 1.0
122    iv % tamdar_sfc_ef_q = 1.0
124    !----------------------------------------------------------------------
125    ! [1.0] Allocate innovation vector and observation structures:
126    !----------------------------------------------------------------------
127    call da_allocate_y(iv, ob)
129    !----------------------------------------------------------------------
130    ! [2.0] Transfer observations:
131    !----------------------------------------------------------------------
133    ! [2.1] Transfer surface obs:
135    if (iv%info(synop)%nlocal > 0) then
136       do n = 1, iv%info(synop)%nlocal
137          ob % synop(n) % u = iv % synop(n) % u % inv
138          ob % synop(n) % v = iv % synop(n) % v % inv
139          ob % synop(n) % t = iv % synop(n) % t % inv
140          ob % synop(n) % q = iv % synop(n) % q % inv
141          ob % synop(n) % p = iv % synop(n) % p % inv
143         if ( q_error_options == 1 ) then
144          ! Calculate q error from rh error:
146          if (.not. present(uvq_direct) .or. (present(uvq_direct) .and. (.not. uvq_direct))) then
147          rh_error = iv%synop(n)%q%error ! q error is rh at this stage!
149          ! if((ob % synop(n) % p > iv%ptop) .AND. &
150          !    (ob % synop(n) % t > 100.0) .AND. &
151          !    (ob % synop(n) % q > 0.0) .AND. &
152          !    (iv % synop(n) % p % qc >= obs_qc_pointer) .and. &
153          !    (iv % synop(n) % t % qc >= obs_qc_pointer) .and. &
154          !    (iv % synop(n) % q % qc >= obs_qc_pointer)) then
155          call da_get_q_error(ob % synop(n) % p, &
156                               ob % synop(n) % t, &
157                               ob % synop(n) % q, &
158                               iv % synop(n) % t % error, &
159                               rh_error, iv % synop(n) % q % error)
160          if (iv%synop(n)% q % error == missing_r) iv%synop(n)% q % qc = missing_data
162          ! end if
163          end if
164         end if !q_error_options
165       end do      
166    end if
168    ! [2.2] Transfer metar obs:
170    if (iv%info(metar)%nlocal > 0) then
171       do n = 1, iv%info(metar)%nlocal
172          ob % metar(n) % u = iv % metar(n) % u % inv
173          ob % metar(n) % v = iv % metar(n) % v % inv
174          ob % metar(n) % t = iv % metar(n) % t % inv
175          ob % metar(n) % q = iv % metar(n) % q % inv
176          ob % metar(n) % p = iv % metar(n) % p % inv
178          ! Calculate q error from rh error:
180          if (.not. present(uvq_direct) .or. (present(uvq_direct) .and. (.not. uvq_direct))) then
181          rh_error = iv%metar(n)%q%error ! q error is rh at this stage!
182          call da_get_q_error(iv % metar(n) % p % inv, &
183                               ob % metar(n) % t, &
184                               ob % metar(n) % q, &
185                               iv % metar(n) % t % error, &
186                               rh_error, q_error)
187          iv % metar(n) % q % error = q_error
188          if (iv%metar(n)% q % error == missing_r) &
189             iv%metar(n)% q % qc = missing_data
190          end if
191       end do
192    end if
194    ! [2.2] Transfer ships obs:
196    if (iv%info(ships)%nlocal > 0) then   
197       do n = 1, iv%info(ships)%nlocal
198          ob % ships(n) % u = iv % ships(n) % u % inv
199          ob % ships(n) % v = iv % ships(n) % v % inv
200          ob % ships(n) % t = iv % ships(n) % t % inv
201          ob % ships(n) % q = iv % ships(n) % q % inv
202          ob % ships(n) % p = iv % ships(n) % p % inv
204          ! Calculate q error from rh error:
206          if (.not. present(uvq_direct) .or. (present(uvq_direct) .and. (.not. uvq_direct))) then
207          rh_error = iv%ships(n)%q%error ! q error is rh at this stage!
208          call da_get_q_error(iv % ships(n) % p % inv, &
209                               ob % ships(n) % t, &
210                               ob % ships(n) % q, &
211                               iv % ships(n) % t % error, &
212                               rh_error, q_error)
213          iv % ships(n) % q % error = q_error
215          if(iv%ships(n)% q % error == missing_r) iv%ships(n)% q % qc = missing_data
216          end if
217       end do
218       
219    end if
221    ! [2.4.1] Transfer Geo. AMVs Obs:
223    if (iv%info(geoamv)%nlocal > 0) then
224       do n = 1, iv%info(geoamv)%nlocal
225          do k = 1, iv%info(geoamv)%levels(n)
226             ob % geoamv(n) % u(k) = iv % geoamv(n) % u(k) % inv
227             ob % geoamv(n) % v(k) = iv % geoamv(n) % v(k) % inv
228          end do
229       end do
230    end if
232    ! [2.4.2] Transfer  Polar AMVs Obs:
234    if (iv%info(polaramv)%nlocal > 0) then
235       do n = 1, iv%info(polaramv)%nlocal
236          do k = 1,iv%info(polaramv)%levels(n)
237             ob % polaramv(n) % u(k) = iv % polaramv(n) % u(k) % inv
238             ob % polaramv(n) % v(k) = iv % polaramv(n) % v(k) % inv
239          end do
240       end do
241    end if
243    ! [2.5] Transfer gpspw obs:
245    if (iv%info(gpspw)%nlocal > 0) then
246       do n = 1, iv%info(gpspw)%nlocal
247          ob % gpspw(n) % tpw = iv % gpspw(n) % tpw % inv
248       end do
250    end if
252    ! [2.6] Transfer GPS REF obs:
254    if (iv%info(gpsref)%nlocal > 0) then
255       do n = 1, iv%info(gpsref)%nlocal
256          do k = 1, iv%info(gpsref)%levels(n)
257 ! ............................................................................
258 ! Convert the geometric height to the geopotential height for GPSREF data
259 ! because GPSRO used geometric height for impact parameter, and the WRF model
260 ! always uses the geopotential height as the vertical coordinate for P, T, q.
261 ! YRG (10/05/2010):
262             geometric_h = iv%gpsref(n)%h(k) / 1000.0
263             call da_msl2geo1 (geometric_h, iv%info(gpsref)%lat(1,n),  &
264                            iv%info(gpsref)%lon(1,n), geopotential_h)
265             iv % gpsref(n) % h(k) = geopotential_h * 1000.0
266 !            write (999,'("n=",i6," k=",i3,2x,"gmh, gph, lat, lon:",4f15.5)') &
267 !                          n, k, geometric_h*1000.0, iv % gpsref(n) % h(k), &
268 !                          iv%info(gpsref)%lat(1,n), iv%info(gpsref)%lon(1,n)
269 ! ............................................................................
270             ob % gpsref(n) % ref(k) = iv % gpsref(n) % ref(k) % inv
271             ob % gpsref(n) %   p(k) = iv % gpsref(n) %   p(k) % inv
272             ob % gpsref(n) %   t(k) = iv % gpsref(n) %   t(k) % inv
273             ob % gpsref(n) %   q(k) = iv % gpsref(n) %   q(k) % inv
274          end do
275       end do
276    end if
278    ! [2.6.1] Transfer GPS EPH obs:
280    if (iv%info(gpseph)%nlocal > 0) then
281       do n = 1, iv%info(gpseph)%nlocal
282          do k = 1, iv%info(gpseph)%levels(n)
283             ob % gpseph(n) % eph(k) = iv % gpseph(n) % eph(k) % inv
284          end do
285       end do
286    end if
288    ! [2.7] Transfer sonde obs:
290    if (iv%info(sound)%nlocal > 0) then
291       do n = 1, iv%info(sound)%nlocal
292          do k = 1, iv%info(sound)%levels(n)
293             ob % sound(n) % u(k) = iv % sound(n) % u(k) % inv
294             ob % sound(n) % v(k) = iv % sound(n) % v(k) % inv
295             ob % sound(n) % t(k) = iv % sound(n) % t(k) % inv
296             ob % sound(n) % q(k) = iv % sound(n) % q(k) % inv
298             ! Calculate q error from rh error:
300          if (.not. present(uvq_direct) .or. (present(uvq_direct) .and. (.not. uvq_direct))) then
301             rh_error = iv%sound(n)%q(k)%error ! q error is rh at this stage!
302             call da_get_q_error(iv % sound(n) % p(k), &
303                                  ob % sound(n) % t(k), &
304                                  ob % sound(n) % q(k), &
305                                  iv % sound(n) % t(k) % error, &
306                                  rh_error, q_error)
308             iv % sound(n) % q(k) % error = q_error
309          if (iv%sound(n)% q(k) % error == missing_r) &
310             iv%sound(n)% q(k) % qc = missing_data
311          end if
312          end do
313       end do
314    end if
316    if (iv%info(sonde_sfc)%nlocal > 0) then
317       do n = 1, iv%info(sonde_sfc)%nlocal
318          ob % sonde_sfc(n) % u = iv % sonde_sfc(n) % u % inv
319          ob % sonde_sfc(n) % v = iv % sonde_sfc(n) % v % inv
320          ob % sonde_sfc(n) % t = iv % sonde_sfc(n) % t % inv
321          ob % sonde_sfc(n) % q = iv % sonde_sfc(n) % q % inv
322          ob % sonde_sfc(n) % p = iv % sonde_sfc(n) % p % inv
324          ! Calculate q error from rh error:
326          if (.not. present(uvq_direct) .or. (present(uvq_direct) .and. (.not. uvq_direct))) then
327          rh_error = iv%sonde_sfc(n)%q%error ! q error is rh at this stage!
328          call da_get_q_error(iv % sonde_sfc(n) % p % inv, &
329                               ob % sonde_sfc(n) % t, &
330                               ob % sonde_sfc(n) % q, &
331                               iv % sonde_sfc(n) % t % error, &
332                               rh_error, iv % sonde_sfc(n) % q % error)
333          if (iv%sonde_sfc(n)% q % error == missing_r) &
334             iv%sonde_sfc(n)% q % qc = missing_data
335          end if
336       end do
337    end if
339    ! [2.8] Transfer airep obs:
341    if (iv%info(airep)%nlocal > 0) then
342       do n = 1, iv%info(airep)%nlocal
343          do k = 1, iv%info(airep)%levels(n)
344             ob % airep(n) % u(k) = iv % airep(n) % u(k) % inv
345             ob % airep(n) % v(k) = iv % airep(n) % v(k) % inv
346             ob % airep(n) % t(k) = iv % airep(n) % t(k) % inv
347             ob % airep(n) % q(k) = iv % airep(n) % q(k) % inv
349             if (.not. present(uvq_direct) .or. (present(uvq_direct) .and. (.not. uvq_direct))) then
350                 rh_error = iv%airep(n)%q(k)%error ! q error is rh at this stage!
351                 call da_get_q_error(iv % airep(n) % p(k), &
352                                     ob % airep(n) % t(k), &
353                                     ob % airep(n) % q(k), &
354                                     iv % airep(n) % t(k) % error, &
355                                     rh_error, q_error)
357                 iv % airep(n) % q(k) % error = q_error
358                 if (iv%airep(n)% q(k) % error == missing_r) &
359                     iv%airep(n)% q(k) % qc = missing_data
360             end if
361          end do
362       end do
363    end if
365    ! [2.9] Transfer pilot obs:
367    if (iv%info(pilot)%nlocal > 0) then
368       do n = 1, iv%info(pilot)%nlocal
369          do k = 1, iv%info(pilot)%levels(n)
370             ob % pilot(n) % u(k) = iv % pilot(n) % u(k) % inv
371             ob % pilot(n) % v(k) = iv % pilot(n) % v(k) % inv
372          end do
373       end do
374    end if
376    ! [2.10] Transfer SSM/I obs:SSMI:
378    if (iv%info(ssmi_rv)%nlocal > 0) then
379       do n = 1, iv%info(ssmi_rv)%nlocal
380          ob % ssmi_rv(n) % speed = iv % ssmi_rv(n) % speed % inv
381          ob % ssmi_rv(n) % tpw   = iv % ssmi_rv(n) % tpw % inv
382       end do
383    end if
385    if (iv%info(ssmi_tb)%nlocal > 0) then
386       do n = 1, iv%info(ssmi_tb)%nlocal
387          ob % ssmi_tb(n) % tb19v = iv % ssmi_tb(n) % tb19v % inv
388          ob % ssmi_tb(n) % tb19h = iv % ssmi_tb(n) % tb19h % inv
389          ob % ssmi_tb(n) % tb22v = iv % ssmi_tb(n) % tb22v % inv
390          ob % ssmi_tb(n) % tb37v = iv % ssmi_tb(n) % tb37v % inv
391          ob % ssmi_tb(n) % tb37h = iv % ssmi_tb(n) % tb37h % inv
392          ob % ssmi_tb(n) % tb85v = iv % ssmi_tb(n) % tb85v % inv
393          ob % ssmi_tb(n) % tb85h = iv % ssmi_tb(n) % tb85h % inv
394       end do
395    end if
397    ! [2.11] Transfer satem obs:
399    if (iv%info(satem)%nlocal > 0) then
400       do n = 1, iv%info(satem)%nlocal
401          do k = 1, iv%info(satem)%levels(n)
402             ob % satem(n) % thickness(k) = iv % satem(n) % thickness(k) % inv
403          end do
404       end do
405    end if
406    
407    ! [2.12] Transfer ssmt1 obs:
409    if (iv%info(ssmt1)%nlocal > 0) then
410       do n = 1, iv%info(ssmt1)%nlocal
411          do k = 1, iv%info(ssmt1)%levels(n)
412             ob % ssmt1(n) % t(k) = iv % ssmt1(n) % t(k) % inv
413          end do
414       end do
416    end if
418    ! [2.13] Transfer ssmt2 obs:
420    if (iv%info(ssmt2)%nlocal > 0) then
421       do n = 1, iv%info(ssmt2)%nlocal
422          do k = 1, iv%info(ssmt2)%levels(n)
423             ob % ssmt2(n) % rh(k) = iv % ssmt2(n) % rh(k) % inv
424          end do
425       end do
426    end if
428    ! [2.15] Transfer scatterometer obs:
430    if (iv%info(qscat)%nlocal > 0) then
431       do n = 1, iv%info(qscat)%nlocal
432          ob % qscat(n) % u = iv % qscat(n) % u % inv
433          ob % qscat(n) % v = iv % qscat(n) % v % inv
434       end do     
435    end if
437    ! [2.16] Transfer profiler obs:
439    if (iv%info(profiler)%nlocal > 0) then
440       do n = 1, iv%info(profiler)%nlocal
441          do k = 1, iv%info(profiler)%levels(n)
442             ob % profiler(n) % u(k) = iv % profiler(n) % u(k) % inv
443             ob % profiler(n) % v(k) = iv % profiler(n) % v(k) % inv
444          end do
445       end do
446    end if
448    ! [2.17] Transfer buoy obs:
450    if (iv%info(buoy)%nlocal > 0) then
451       do n = 1, iv%info(buoy)%nlocal
452          ob % buoy(n) % p = iv % buoy(n) % p % inv
453       end do
454       do n = 1, iv%info(buoy)%nlocal
455          ob % buoy(n) % u = iv % buoy(n) % u % inv
456          ob % buoy(n) % v = iv % buoy(n) % v % inv
457          ob % buoy(n) % t = iv % buoy(n) % t % inv
458          ob % buoy(n) % q = iv % buoy(n) % q % inv
460          ! Calculate q error from rh error:
462          if (.not. present(uvq_direct) .or. (present(uvq_direct) .and. (.not. uvq_direct))) then
463          rh_error = iv%buoy(n)%q%error ! q error is rh at this stage!
464          call da_get_q_error(iv % buoy(n) % p % inv, &
465                               ob % buoy(n) % t, &
466                               ob % buoy(n) % q, &
467                               iv % buoy(n) % t % error, &
468                               rh_error, q_error)
469          iv % buoy(n) % q % error = q_error
471          if(iv%buoy (n)% q % error == missing_r) iv%buoy (n)% q % qc = missing_data
472          end if
473       end do
474    end if
476    ! [2.18] Transfer radar obs: this section has been moved to da_fill_obs_structures_radar.inc
478 !   if (iv%info(radar)%nlocal > 0) then
479 !      do n = 1, iv%info(radar)%nlocal
480 !         do k = 1, iv%info(radar)%levels(n)
481 !            ! Copy observation variables:
482 !            ob % radar(n) % rv(k) = iv % radar(n) % rv(k) % inv
483 !           ob % radar(n) % rf(k) = iv % radar(n) % rf(k) % inv
484 !         end do
485 !      end do
486 !   end if
488    ! [2.19] Transfer TC bogus:
490    if (iv%info(bogus)%nlocal > 0) then
491       do n = 1, iv%info(bogus)%nlocal
492          do k = 1, iv%info(bogus)%levels(n)
494             ! Copy observation variables:
496             ob % bogus(n) % u(k) = iv % bogus(n) % u(k) % inv
497             ob % bogus(n) % v(k) = iv % bogus(n) % v(k) % inv
498             ob % bogus(n) % t(k) = iv % bogus(n) % t(k) % inv
499             ob % bogus(n) % q(k) = iv % bogus(n) % q(k) % inv
501             ! Calculate q error from rh error:
502             ! the conversion is not needed for LSAC
503             if (iv%info(bogus)%name(n)(1:4) /= "LSAC") then
504                rh_error = iv%bogus(n)%q(k)%error ! q error is rh at this stage!
505                call da_get_q_error(iv % bogus(n) % p(k), &
506                                    ob % bogus(n) % t(k), &
507                                    ob % bogus(n) % q(k), &
508                                    iv % bogus(n) % t(k) % error, &
509                                    rh_error, q_error)
511                iv % bogus(n) % q(k) % error = q_error
512                if (iv%bogus(n)% q(k) % error == missing_r) &
513                   iv%bogus(n)% q(k) % qc = missing_data
514             end if !LSAC
515          end do
516          ob % bogus(n) % slp = iv % bogus(n) % slp % inv
517       end do
518    end if
520    ! [2.20] Transfer rain obs:
522    if (iv%info(rain)%nlocal > 0) then
523       do n = 1, iv%info(rain)%nlocal
524             ob % rain(n) % rain = iv % rain(n) % rain % inv
525       end do
526    end if
528    ! Transfer AIRS  retrievals:
530    if (iv%info(airsr)%nlocal > 0) then
531       do n = 1, iv%info(airsr)%nlocal
532          do k = 1, iv%info(airsr)%levels(n)
534             ! Copy observation variables:
536             ob % airsr(n) % t(k) = iv % airsr(n) % t(k) % inv
537             ob % airsr(n) % q(k) = iv % airsr(n) % q(k) % inv
539             ! Calculate q error from rh error:
541          if (.not. present(uvq_direct) .or. (present(uvq_direct) .and. (.not. uvq_direct))) then
542             rh_error = iv%airsr(n)%q(k)%error ! q error is rh at this stage!
543             call da_get_q_error(iv % airsr(n) % p(k), &
544                                  ob % airsr(n) % t(k), &
545                                  ob % airsr(n) % q(k), &
546                                  iv % airsr(n) % t(k) % error, &
547                                  rh_error, q_error)
549             iv % airsr(n) % q(k) % error = q_error
550             if (iv%airsr(n)% q(k) % error == missing_r) &
551                iv%airsr(n)% q(k) % qc = missing_data
552          end if
553          end do
554       end do
555    end if
556    if (iv%info(mtgirs)%nlocal > 0) then
557       do n = 1, iv%info(mtgirs)%nlocal
558          do k = 1, iv%info(mtgirs)%levels(n)
559             ob % mtgirs(n) % u(k) = iv % mtgirs(n) % u(k) % inv
560             ob % mtgirs(n) % v(k) = iv % mtgirs(n) % v(k) % inv
561             ob % mtgirs(n) % t(k) = iv % mtgirs(n) % t(k) % inv
562             ob % mtgirs(n) % q(k) = iv % mtgirs(n) % q(k) % inv
563          if (iv%mtgirs(n)% q(k) % error == missing_r) &
564             iv%mtgirs(n)% q(k) % qc = missing_data
565          end do
566       end do
567    end if
569    if (iv%info(tamdar)%nlocal > 0) then
570       do n = 1, iv%info(tamdar)%nlocal
571          do k = 1, iv%info(tamdar)%levels(n)
572             ob % tamdar(n) % u(k) = iv % tamdar(n) % u(k) % inv
573             ob % tamdar(n) % v(k) = iv % tamdar(n) % v(k) % inv
574             ob % tamdar(n) % t(k) = iv % tamdar(n) % t(k) % inv
575             ob % tamdar(n) % q(k) = iv % tamdar(n) % q(k) % inv
577          if (iv%tamdar(n)% u(k) % error == missing_r) &
578             iv%tamdar(n)% u(k) % qc = missing_data
579          if (iv%tamdar(n)% v(k) % error == missing_r) &
580             iv%tamdar(n)% v(k) % qc = missing_data
581          if (iv%tamdar(n)% t(k) % error == missing_r) &
582             iv%tamdar(n)% t(k) % qc = missing_data
584             ! Calculate q error from rh error:
586             rh_error = iv%tamdar(n)%q(k)%error ! q error is rh at this stage!
587             call da_get_q_error(iv % tamdar(n) % p(k), &
588                                 ob % tamdar(n) % t(k), &
589                                 ob % tamdar(n) % q(k), &
590                                 iv % tamdar(n) % t(k) % error, &
591                                 rh_error, q_error)
593             iv % tamdar(n) % q(k) % error = q_error
595          if (iv%tamdar(n)% q(k) % error == missing_r) &
596             iv%tamdar(n)% q(k) % qc = missing_data
597          end do
598       end do
599    end if
601    if (iv%info(tamdar_sfc)%nlocal > 0) then
602       do n = 1, iv%info(tamdar_sfc)%nlocal
604          ob % tamdar_sfc(n) % u = iv % tamdar_sfc(n) % u % inv
605          ob % tamdar_sfc(n) % v = iv % tamdar_sfc(n) % v % inv
606          ob % tamdar_sfc(n) % t = iv % tamdar_sfc(n) % t % inv
607          ob % tamdar_sfc(n) % q = iv % tamdar_sfc(n) % q % inv
608          ob % tamdar_sfc(n) % p = iv % tamdar_sfc(n) % p % inv
610          if (iv%tamdar_sfc(n)% u % error == missing_r) &
611             iv%tamdar_sfc(n)% u % qc = missing_data
612          if (iv%tamdar_sfc(n)% v % error == missing_r) &
613             iv%tamdar_sfc(n)% v % qc = missing_data
614          if (iv%tamdar_sfc(n)% t % error == missing_r) &
615             iv%tamdar_sfc(n)% t % qc = missing_data
617          ! Calculate q error from rh error:
619          rh_error = iv%tamdar_sfc(n)%q%error ! q error is rh at this stage!
620          call da_get_q_error(iv % tamdar_sfc(n) % p % inv, &
621                               ob % tamdar_sfc(n) % t, &
622                               ob % tamdar_sfc(n) % q, &
623                               iv % tamdar_sfc(n) % t % error, &
624                               rh_error, iv % tamdar_sfc(n) % q % error)
625          if (iv%tamdar_sfc(n)% q % error == missing_r) &
626             iv%tamdar_sfc(n)% q % qc = missing_data
627       end do
628    end if
630    if (trace_use) call da_trace_exit("da_fill_obs_structures")
632 end subroutine da_fill_obs_structures