da_fill_obs_structures.inc

References to this file elsewhere.
1 subroutine da_fill_obs_structures(iv, ob)
2 
3    !----------------------------------------------------------------------------   
4    ! Purpose: Allocates observation structure and fills it from iv.
5    !----------------------------------------------------------------------------   
6 
7    implicit none
8 
9    type (iv_type), intent(inout) :: iv   ! Obs and header structure.
10    type (y_type), intent(out)    :: ob   ! (Smaller) observation structure.
11 
12    integer :: n, k     ! Loop counters.
13    real    :: rh_error ! RH obs. error.
14    real    :: q_error  ! q obs. error.
15    integer :: i,j
16    logical :: outside
17 
18    if (trace_use) call da_trace_entry("da_fill_obs_structures")
19 
20    !---------------------------------------------------------------------------
21    ! Initialise obs error factors (which will be overwritten in use_obs_errfac)
22    !---------------------------------------------------------------------------
23 
24    iv % synop_ef_u = 1.0
25    iv % synop_ef_v = 1.0
26    iv % synop_ef_t = 1.0
27    iv % synop_ef_p = 1.0
28    iv % synop_ef_q = 1.0
29 
30    iv % metar_ef_u = 1.0
31    iv % metar_ef_v = 1.0
32    iv % metar_ef_t = 1.0
33    iv % metar_ef_p = 1.0
34    iv % metar_ef_q = 1.0
35 
36    iv % ships_ef_u = 1.0
37    iv % ships_ef_v = 1.0
38    iv % ships_ef_t = 1.0
39    iv % ships_ef_p = 1.0
40    iv % ships_ef_q = 1.0
41 
42    iv % geoamv_ef_u = 1.0
43    iv % geoamv_ef_v = 1.0
44 
45    iv % polaramv_ef_u = 1.0
46    iv % polaramv_ef_v = 1.0
47 
48    iv % gpspw_ef_tpw = 1.0
49 
50    iv % gpsref_ef_ref = 1.0
51    iv % gpsref_ef_p = 1.0
52    iv % gpsref_ef_t = 1.0
53    iv % gpsref_ef_q = 1.0
54 
55    iv % sound_ef_u = 1.0
56    iv % sound_ef_v = 1.0
57    iv % sound_ef_t = 1.0
58    iv % sound_ef_q = 1.0
59 
60    iv % airep_ef_u = 1.0
61    iv % airep_ef_v = 1.0
62    iv % airep_ef_t = 1.0
63 
64    iv % pilot_ef_u = 1.0
65    iv % pilot_ef_v = 1.0
66 
67    iv % ssmir_ef_speed = 1.0
68    iv % ssmir_ef_tpw = 1.0
69 
70    iv % satem_ef_thickness = 1.0
71 
72    iv % ssmt1_ef_t = 1.0
73 
74    iv % ssmt2_ef_rh = 1.0
75 
76    iv % qscat_ef_u = 1.0
77    iv % qscat_ef_v = 1.0
78 
79    iv % profiler_ef_u = 1.0
80    iv % profiler_ef_v = 1.0
81    
82    iv % buoy_ef_u = 1.0
83    iv % buoy_ef_v = 1.0
84    iv % buoy_ef_t = 1.0
85    iv % buoy_ef_p = 1.0
86    iv % buoy_ef_q = 1.0
87 
88    iv % radar_ef_rv = 1.0
89    iv % radar_ef_rf = 1.0
90 
91    iv % bogus_ef_u = 1.0
92    iv % bogus_ef_v = 1.0
93    iv % bogus_ef_t = 1.0
94    iv % bogus_ef_p = 1.0
95    iv % bogus_ef_q = 1.0
96    iv % bogus_ef_slp = 1.0
97 
98    iv % airsr_ef_t = 1.0
99    iv % airsr_ef_q = 1.0
100 
101    !-------------------------------
102    ! Find out if pseudo ob is local
103    !-------------------------------
104 
105    if (num_pseudo > 0) then
106       outside = .false.
107       i = int(pseudo_x)
108       j = int(pseudo_y)
109 
110       if (fg_format == fg_format_kma_global) then
111          if ((j < jts-1) .or. (j > jte)) outside = .true.
112       else
113          if ((i < ids)   .or. (i >= ide) .or. (j < jds)   .or. (j >= jde)) outside = .true.
114          if ((i < its-1) .or. (i >  ite) .or. (j < jts-1) .or. (j >  jte)) outside = .true.
115       end if
116    else
117       outside = .true.
118    end if
119 
120    if (outside) then
121       iv%info(pseudo)%nlocal          = 0
122       ob%nlocal(pseudo)               = 0
123       iv%info(pseudo)%ntotal          = 0
124       iv%info(pseudo)%plocal(iv%time) = 0
125    else
126       iv%info(pseudo)%nlocal          = num_pseudo
127       ob%nlocal(pseudo)               = num_pseudo
128       iv%info(pseudo)%ntotal          = num_pseudo
129       iv%info(pseudo)%plocal(iv%time) = num_pseudo
130    end if
131 
132    !---------------------------------------------------------------------------
133    ! [1.0] Allocate innovation vector and observation structures:
134    !---------------------------------------------------------------------------
135    call da_allocate_y(iv, ob)
136 
137    !----------------------------------------------------------------------
138    ! [2.0] Transfer observations:
139    !----------------------------------------------------------------------
140 
141    ! [2.1] Transfer surface obs:
142 
143    if (iv%info(synop)%nlocal > 0) then
144       do n = 1, iv%info(synop)%nlocal
145          ob % synop(n) % u = iv % synop(n) % u % inv
146          ob % synop(n) % v = iv % synop(n) % v % inv
147          ob % synop(n) % t = iv % synop(n) % t % inv
148          ob % synop(n) % q = iv % synop(n) % q % inv
149          ob % synop(n) % p = iv % synop(n) % p % inv
150 
151          ! Calculate q error from rh error:
152 
153          rh_error = iv%synop(n)%q%error ! q error is rh at this stage!
154 
155          ! if((ob % synop(n) % p > iv%ptop) .AND. &
156          !    (ob % synop(n) % t > 100.0) .AND. &
157          !    (ob % synop(n) % q > 0.0) .AND. &
158          !    (iv % synop(n) % p % qc >= obs_qc_pointer) .and. &
159          !    (iv % synop(n) % t % qc >= obs_qc_pointer) .and. &
160          !    (iv % synop(n) % q % qc >= obs_qc_pointer)) then
161          call da_get_q_error(ob % synop(n) % p, &
162                               ob % synop(n) % t, &
163                               ob % synop(n) % q, &
164                               iv % synop(n) % t % error, &
165                               rh_error, iv % synop(n) % q % error)
166          if (iv%synop(n)% q % error == missing_r) iv%synop(n)% q % qc = missing_data
167 
168          ! end if
169       end do      
170    end if
171 
172    ! [2.2] Transfer metar obs:
173 
174    if (iv%info(metar)%nlocal > 0) then
175       do n = 1, iv%info(metar)%nlocal
176          ob % metar(n) % u = iv % metar(n) % u % inv
177          ob % metar(n) % v = iv % metar(n) % v % inv
178          ob % metar(n) % t = iv % metar(n) % t % inv
179          ob % metar(n) % q = iv % metar(n) % q % inv
180          ob % metar(n) % p = iv % metar(n) % p % inv
181 
182          ! Calculate q error from rh error:
183 
184          rh_error = iv%metar(n)%q%error ! q error is rh at this stage!
185          call da_get_q_error(iv % metar(n) % p % inv, &
186                               ob % metar(n) % t, &
187                               ob % metar(n) % q, &
188                               iv % metar(n) % t % error, &
189                               rh_error, q_error)
190          iv % metar(n) % q % error = q_error
191          if (iv%metar(n)% q % error == missing_r) &
192             iv%metar(n)% q % qc = missing_data
193       end do
194    end if
195 
196    ! [2.2] Transfer ships obs:
197 
198    if (iv%info(ships)%nlocal > 0) then   
199       do n = 1, iv%info(ships)%nlocal
200          ob % ships(n) % u = iv % ships(n) % u % inv
201          ob % ships(n) % v = iv % ships(n) % v % inv
202          ob % ships(n) % t = iv % ships(n) % t % inv
203          ob % ships(n) % q = iv % ships(n) % q % inv
204          ob % ships(n) % p = iv % ships(n) % p % inv
205 
206          ! Calculate q error from rh error:
207 
208          rh_error = iv%ships(n)%q%error ! q error is rh at this stage!
209          call da_get_q_error(iv % ships(n) % p % inv, &
210                               ob % ships(n) % t, &
211                               ob % ships(n) % q, &
212                               iv % ships(n) % t % error, &
213                               rh_error, q_error)
214          iv % ships(n) % q % error = q_error
215 
216          if(iv%ships(n)% q % error == missing_r) iv%ships(n)% q % qc = missing_data
217       end do
218       
219    end if
220 
221    ! [2.4.1] Transfer Geo. AMVs Obs:
222 
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
231 
232    ! [2.4.2] Transfer  Polar AMVs Obs:
233 
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
242 
243    ! [2.5] Transfer gpspw obs:
244 
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
249 
250    end if
251 
252    ! [2.6] Transfer GPS REF obs:
253 
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             ob % gpsref(n) % ref(k) = iv % gpsref(n) % ref(k) % inv
258             ob % gpsref(n) %   p(k) = iv % gpsref(n) %   p(k) % inv
259             ob % gpsref(n) %   t(k) = iv % gpsref(n) %   t(k) % inv
260             ob % gpsref(n) %   q(k) = iv % gpsref(n) %   q(k) % inv
261          end do
262       end do
263    end if
264 
265    ! [2.7] Transfer sonde obs:
266 
267    if (iv%info(sound)%nlocal > 0) then
268       do n = 1, iv%info(sound)%nlocal
269          do k = 1, iv%info(sound)%levels(n)
270             ob % sound(n) % u(k) = iv % sound(n) % u(k) % inv
271             ob % sound(n) % v(k) = iv % sound(n) % v(k) % inv
272             ob % sound(n) % t(k) = iv % sound(n) % t(k) % inv
273             ob % sound(n) % q(k) = iv % sound(n) % q(k) % inv
274 
275             ! Calculate q error from rh error:
276 
277             rh_error = iv%sound(n)%q(k)%error ! q error is rh at this stage!
278             call da_get_q_error(iv % sound(n) % p(k), &
279                                  ob % sound(n) % t(k), &
280                                  ob % sound(n) % q(k), &
281                                  iv % sound(n) % t(k) % error, &
282                                  rh_error, q_error)
283 
284             iv % sound(n) % q(k) % error = q_error
285          if (iv%sound(n)% q(k) % error == missing_r) &
286             iv%sound(n)% q(k) % qc = missing_data
287          end do
288          ob % sonde_sfc(n) % u = iv % sonde_sfc(n) % u % inv
289          ob % sonde_sfc(n) % v = iv % sonde_sfc(n) % v % inv
290          ob % sonde_sfc(n) % t = iv % sonde_sfc(n) % t % inv
291          ob % sonde_sfc(n) % q = iv % sonde_sfc(n) % q % inv
292          ob % sonde_sfc(n) % p = iv % sonde_sfc(n) % p % inv
293 
294          ! Calculate q error from rh error:
295 
296          rh_error = iv%sonde_sfc(n)%q%error ! q error is rh at this stage!
297          call da_get_q_error(iv % sonde_sfc(n) % p % inv, &
298                               ob % sonde_sfc(n) % t, &
299                               ob % sonde_sfc(n) % q, &
300                               iv % sonde_sfc(n) % t % error, &
301                               rh_error, iv % sonde_sfc(n) % q % error)
302          if (iv%sonde_sfc(n)% q % error == missing_r) &
303             iv%sonde_sfc(n)% q % qc = missing_data
304       end do
305    end if
306 
307    ! [2.8] Transfer airep obs:
308 
309    if (iv%info(airep)%nlocal > 0) then
310       do n = 1, iv%info(airep)%nlocal
311          do k = 1, iv%info(airep)%levels(n)
312             ob % airep(n) % u(k) = iv % airep(n) % u(k) % inv
313             ob % airep(n) % v(k) = iv % airep(n) % v(k) % inv
314             ob % airep(n) % t(k) = iv % airep(n) % t(k) % inv
315          end do
316       end do
317    end if
318 
319    ! [2.9] Transfer pilot obs:
320 
321    if (iv%info(pilot)%nlocal > 0) then
322       do n = 1, iv%info(pilot)%nlocal
323          do k = 1, iv%info(pilot)%levels(n)
324             ob % pilot(n) % u(k) = iv % pilot(n) % u(k) % inv
325             ob % pilot(n) % v(k) = iv % pilot(n) % v(k) % inv
326          end do
327       end do
328    end if
329 
330    ! [2.10] Transfer SSM/I obs:SSMI:
331 
332    if (iv%info(ssmi_rv)%nlocal > 0) then
333       do n = 1, iv%info(ssmi_rv)%nlocal
334          ob % ssmi_rv(n) % speed = iv % ssmi_rv(n) % speed % inv
335          ob % ssmi_rv(n) % tpw   = iv % ssmi_rv(n) % tpw % inv
336       end do
337    end if
338 
339    if (iv%info(ssmi_tb)%nlocal > 0) then
340       do n = 1, iv%info(ssmi_tb)%nlocal
341          ob % ssmi_tb(n) % tb19v = iv % ssmi_tb(n) % tb19v % inv
342          ob % ssmi_tb(n) % tb19h = iv % ssmi_tb(n) % tb19h % inv
343          ob % ssmi_tb(n) % tb22v = iv % ssmi_tb(n) % tb22v % inv
344          ob % ssmi_tb(n) % tb37v = iv % ssmi_tb(n) % tb37v % inv
345          ob % ssmi_tb(n) % tb37h = iv % ssmi_tb(n) % tb37h % inv
346          ob % ssmi_tb(n) % tb85v = iv % ssmi_tb(n) % tb85v % inv
347          ob % ssmi_tb(n) % tb85h = iv % ssmi_tb(n) % tb85h % inv
348       end do
349    end if
350 
351    ! [2.11] Transfer satem obs:
352 
353    if (iv%info(satem)%nlocal > 0) then
354       do n = 1, iv%info(satem)%nlocal
355          do k = 1, iv%info(satem)%levels(n)
356             ob % satem(n) % thickness(k) = iv % satem(n) % thickness(k) % inv
357          end do
358       end do
359    end if
360    
361    ! [2.12] Transfer ssmt1 obs:
362 
363    if (iv%info(ssmt1)%nlocal > 0) then
364       do n = 1, iv%info(ssmt1)%nlocal
365          do k = 1, iv%info(ssmt1)%levels(n)
366             ob % ssmt1(n) % t(k) = iv % ssmt1(n) % t(k) % inv
367          end do
368       end do
369 
370    end if
371 
372    ! [2.13] Transfer ssmt2 obs:
373 
374    if (iv%info(ssmt2)%nlocal > 0) then
375       do n = 1, iv%info(ssmt2)%nlocal
376          do k = 1, iv%info(ssmt2)%levels(n)
377             ob % ssmt2(n) % rh(k) = iv % ssmt2(n) % rh(k) % inv
378          end do
379       end do
380    end if
381    
382    ! [2.14] Setup pseudo observations:
383 
384    if (iv%info(pseudo)%nlocal > 0) call da_setup_pseudo_obs(iv, ob)
385 
386    ! [2.15] Transfer scatterometer obs:
387 
388    if (iv%info(qscat)%nlocal > 0) then
389       do n = 1, iv%info(qscat)%nlocal
390          ob % qscat(n) % u = iv % qscat(n) % u % inv
391          ob % qscat(n) % v = iv % qscat(n) % v % inv
392       end do     
393    end if
394 
395    ! [2.16] Transfer profiler obs:
396 
397    if (iv%info(profiler)%nlocal > 0) then
398       do n = 1, iv%info(profiler)%nlocal
399          do k = 1, iv%info(profiler)%levels(n)
400             ob % profiler(n) % u(k) = iv % profiler(n) % u(k) % inv
401             ob % profiler(n) % v(k) = iv % profiler(n) % v(k) % inv
402          end do
403       end do
404    end if
405 
406    ! [2.17] Transfer buoy obs:
407 
408    if (iv%info(buoy)%nlocal > 0) then
409       do n = 1, iv%info(buoy)%nlocal
410          ob % buoy(n) % p = iv % buoy(n) % p % inv
411       end do
412       do n = 1, iv%info(buoy)%nlocal
413          ob % buoy(n) % u = iv % buoy(n) % u % inv
414          ob % buoy(n) % v = iv % buoy(n) % v % inv
415          ob % buoy(n) % t = iv % buoy(n) % t % inv
416          ob % buoy(n) % q = iv % buoy(n) % q % inv
417 
418          ! Calculate q error from rh error:
419 
420          rh_error = iv%buoy(n)%q%error ! q error is rh at this stage!
421          call da_get_q_error(iv % buoy(n) % p % inv, &
422                               ob % buoy(n) % t, &
423                               ob % buoy(n) % q, &
424                               iv % buoy(n) % t % error, &
425                               rh_error, q_error)
426          iv % buoy(n) % q % error = q_error
427 
428          if(iv%buoy (n)% q % error == missing_r) iv%buoy (n)% q % qc = missing_data
429       end do
430    end if
431 
432    ! [2.18] Transfer radar obs:
433 
434    if (iv%info(radar)%nlocal > 0) then
435       do n = 1, iv%info(radar)%nlocal
436          do k = 1, iv%info(radar)%levels(n)
437             ! Copy observation variables:
438             ob % radar(n) % rv(k) = iv % radar(n) % rv(k) % inv
439            ob % radar(n) % rf(k) = iv % radar(n) % rf(k) % inv
440          end do
441       end do
442    end if
443 
444    ! [2.19] Transfer TC bogus:
445 
446    if (iv%info(bogus)%nlocal > 0) then
447       do n = 1, iv%info(bogus)%nlocal
448          do k = 1, iv%info(bogus)%levels(n)
449 
450             ! Copy observation variables:
451 
452             ob % bogus(n) % u(k) = iv % bogus(n) % u(k) % inv
453             ob % bogus(n) % v(k) = iv % bogus(n) % v(k) % inv
454             ob % bogus(n) % t(k) = iv % bogus(n) % t(k) % inv
455             ob % bogus(n) % q(k) = iv % bogus(n) % q(k) % inv
456 
457             ! Calculate q error from rh error:
458 
459             rh_error = iv%bogus(n)%q(k)%error ! q error is rh at this stage!
460             call da_get_q_error(iv % bogus(n) % p(k), &
461                                  ob % bogus(n) % t(k), &
462                                  ob % bogus(n) % q(k), &
463                                  iv % bogus(n) % t(k) % error, &
464                                  rh_error, q_error)
465 
466             iv % bogus(n) % q(k) % error = q_error
467             if (iv%bogus(n)% q(k) % error == missing_r) &
468                iv%bogus(n)% q(k) % qc = missing_data
469          end do
470          ob % bogus(n) % slp = iv % bogus(n) % slp % inv
471       end do
472    end if
473 
474    ! Transfer AIRS  retrievals:
475 
476    if (iv%info(airsr)%nlocal > 0) then
477       do n = 1, iv%info(airsr)%nlocal
478          do k = 1, iv%info(airsr)%levels(n)
479 
480             ! Copy observation variables:
481 
482             ob % airsr(n) % t(k) = iv % airsr(n) % t(k) % inv
483             ob % airsr(n) % q(k) = iv % airsr(n) % q(k) % inv
484 
485             ! Calculate q error from rh error:
486 
487             rh_error = iv%airsr(n)%q(k)%error ! q error is rh at this stage!
488             call da_get_q_error(iv % airsr(n) % p(k), &
489                                  ob % airsr(n) % t(k), &
490                                  ob % airsr(n) % q(k), &
491                                  iv % airsr(n) % t(k) % error, &
492                                  rh_error, q_error)
493 
494             iv % airsr(n) % q(k) % error = q_error
495             if (iv%airsr(n)% q(k) % error == missing_r) &
496                iv%airsr(n)% q(k) % qc = missing_data
497          end do
498       end do
499    end if
500 
501    if (trace_use) call da_trace_exit("da_fill_obs_structures")
502 
503 end subroutine da_fill_obs_structures
504 
505