Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / var / da / da_obs / da_use_obs_errfac.inc
blob9842520bd262d6b82a5881f9a4041264427db75b
1 subroutine da_use_obs_errfac(iv)
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.
11    integer                       :: n, k            ! Loop counters.
12    real                          :: d1, d2, d3, d4  ! Dummy values.
14    if (trace_use) call da_trace_entry("da_use_obs_errfac")
16    !----------------------------------------------------------------------
17    ! [2.0] Scale observation errors:
18    !-------------------------------------------------------------------
20    ! [2.1] Transfer surface obs:
22    call da_read_errfac('synop', iv % synop_ef_u, &
23                         iv % synop_ef_v, iv % synop_ef_t, &
24                         iv % synop_ef_p, iv % synop_ef_q)                           
26    if (iv%info(synop)%nlocal > 0) then
27       do n = 1, iv%info(synop)%nlocal
28          iv % synop(n) % u % error = iv % synop(n) % u % error * iv % synop_ef_u
29          iv % synop(n) % v % error = iv % synop(n) % v % error * iv % synop_ef_v
30          iv % synop(n) % t % error = iv % synop(n) % t % error * iv % synop_ef_t
31          iv % synop(n) % p % error = iv % synop(n) % p % error * iv % synop_ef_p
32          iv % synop(n) % q % error = iv % synop(n) % q % error * iv % synop_ef_q
33       end do
34    end if
36    ! [2.2] Transfer metar obs:
39    call da_read_errfac('metar', iv % metar_ef_u, &
40                         iv % metar_ef_v, iv % metar_ef_t, &
41                         iv % metar_ef_p, iv % metar_ef_q)
42                            
43    if (iv%info(metar)%nlocal > 0) then
44       do n = 1, iv%info(metar)%nlocal
45          iv % metar(n) % u % error = iv % metar(n) % u % error * iv % metar_ef_u
46          iv % metar(n) % v % error = iv % metar(n) % v % error * iv % metar_ef_v
47          iv % metar(n) % t % error = iv % metar(n) % t % error * iv % metar_ef_t
48          iv % metar(n) % p % error = iv % metar(n) % p % error * iv % metar_ef_p
49          iv % metar(n) % q % error = iv % metar(n) % q % error * iv % metar_ef_q
50       end do
51    end if
53    ! [2.2] Transfer ships obs:
55       
56    call da_read_errfac('ships', iv % ships_ef_u, &
57                         iv % ships_ef_v, iv % ships_ef_t, &
58                         iv % ships_ef_p, iv % ships_ef_q)
59                            
60    if (iv%info(ships)%nlocal > 0) then
61       do n = 1, iv%info(ships)%nlocal
62          iv % ships(n) % u % error = iv % ships(n) % u % error * iv % ships_ef_u
63          iv % ships(n) % v % error = iv % ships(n) % v % error * iv % ships_ef_v
64          iv % ships(n) % t % error = iv % ships(n) % t % error * iv % ships_ef_t
65          iv % ships(n) % p % error = iv % ships(n) % p % error * iv % ships_ef_p
66          iv % ships(n) % q % error = iv % ships(n) % q % error * iv % ships_ef_q
67       end do
68    end if
70    ! [2.4.1] Transfer Geo. AMVs Obs:
71    
72    
73    call da_read_errfac('geoamv', iv % geoamv_ef_u, iv % geoamv_ef_v, d1, d2, d3)
75    if (iv%info(geoamv)%nlocal > 0) then
76       do n = 1, iv%info(geoamv)%nlocal
77         do k = 1, iv%info(geoamv)%levels(n)
78          iv % geoamv(n) % u(k) % error = iv % geoamv(n) % u(k) % error * iv % geoamv_ef_u
79          iv % geoamv(n) % v(k) % error = iv % geoamv(n) % v(k) % error * iv % geoamv_ef_v
80         end do
81       end do 
82    end if
84    ! [2.4.2] Transfer Polar AMVs Obs:
87    call da_read_errfac('polaramv', iv % polaramv_ef_u, iv % polaramv_ef_v, d1, d2, d3)
89    if (iv%info(polaramv)%nlocal > 0) then 
90       do n = 1, iv%info(polaramv)%nlocal
91         do k = 1, iv%info(polaramv)%levels(n)
92          iv % polaramv(n) % u(k) % error = iv % polaramv(n) % u(k) % error * iv % polaramv_ef_u
93          iv % polaramv(n) % v(k) % error = iv % polaramv(n) % v(k) % error * iv % polaramv_ef_v
94         end do
95       end do
96    end if
99    ! [2.5] Transfer gpspw obs:
102    call da_read_errfac('gpspw', iv % gpspw_ef_tpw, d1, d2, d3, d4)
104    if (iv%info(gpspw)%nlocal > 0) then
105       do n = 1, iv%info(gpspw)%nlocal
106          iv % gpspw(n) % tpw % error = iv % gpspw(n) % tpw % error * &
107                                        iv % gpspw_ef_tpw
109       end do
110    end if
112 ! [2.5.1] Transfer gpsref obs:
114    call da_read_errfac('gpsre', iv % gpsref_ef_ref, d1, d2, d3, d4)
116    if (iv%info(gpsref)%nlocal > 0) then
117       do n = 1, iv%info(gpsref)%nlocal
118          do k = 1, iv%info(gpsref)%levels(n)
119             iv % gpsref(n) % ref(k) % error = iv % gpsref(n) % ref(k) % error * &
120                                               iv % gpsref_ef_ref
121          enddo
122       end do
123    end if
125 ! [2.5.2] Transfer gpseph obs:
127    call da_read_errfac('gpsep', iv % gpseph_ef_eph, d1, d2, d3, d4)
129    if (iv%info(gpseph)%nlocal > 0) then
130       do n = 1, iv%info(gpseph)%nlocal
131          do k = 1, iv%info(gpseph)%levels(n)
132             if ( iv % gpseph(n) % eph(k) % error > 0.0 ) then
133                !only apply the error factor when the error is not missing
134                iv % gpseph(n) % eph(k) % error = iv % gpseph(n) % eph(k) % error * &
135                                                  iv % gpseph_ef_eph
136             end if
137          enddo
138       end do
139    end if
141    ! [2.6] Transfer sonde obs:
144    call da_read_errfac('sound', iv % sound_ef_u, iv % sound_ef_v, &
145                         iv % sound_ef_t, iv % sound_ef_q, d1)
147    if (iv%info(sound)%nlocal > 0) then
148       do n = 1, iv%info(sound)%nlocal
149          do k = 1, iv%info(sound)%levels(n)
150             iv % sound(n) % u(k) % error = iv % sound(n) % u(k) % error * &
151                                            iv % sound_ef_u
152             iv % sound(n) % v(k) % error = iv % sound(n) % v(k) % error * &
153                                            iv % sound_ef_v
154             iv % sound(n) % t(k) % error = iv % sound(n) % t(k) % error * &
155                                            iv % sound_ef_t
156             iv % sound(n) % q(k) % error = iv % sound(n) % q(k) % error * &
157                                            iv % sound_ef_q
158          end do
159       end do
160    end if
162    if (iv%info(sonde_sfc)%nlocal > 0) then
163       do n = 1, iv%info(sonde_sfc)%nlocal
164          iv % sonde_sfc(n) % u % error = iv % sonde_sfc(n) % u % error * iv % synop_ef_u
165          iv % sonde_sfc(n) % v % error = iv % sonde_sfc(n) % v % error * iv % synop_ef_v
166          iv % sonde_sfc(n) % t % error = iv % sonde_sfc(n) % t % error * iv % synop_ef_t
167          iv % sonde_sfc(n) % p % error = iv % sonde_sfc(n) % p % error * iv % synop_ef_p
168          iv % sonde_sfc(n) % q % error = iv % sonde_sfc(n) % q % error * iv % synop_ef_q
169       end do
170    end if
172    ! [2.7] Transfer airep obs:
174    
175    call da_read_errfac('airep', iv % airep_ef_u, iv % airep_ef_v, &
176                         iv % airep_ef_t, iv % airep_ef_q, d1)
178    if (iv%info(airep)%nlocal > 0) then
179       do n = 1, iv%info(airep)%nlocal
180          do k = 1, iv%info(airep)%levels(n)
181             iv % airep(n) % u(k) % error = iv % airep(n) % u(k) % error * &
182                                            iv % airep_ef_u
183             iv % airep(n) % v(k) % error = iv % airep(n) % v(k) % error * &
184                                            iv % airep_ef_v
185             iv % airep(n) % t(k) % error = iv % airep(n) % t(k) % error * &
186                                            iv % airep_ef_t
187             iv % airep(n) % q(k) % error = iv % airep(n) % q(k) % error * &
188                                            iv % airep_ef_q
189          end do
190       end do
191    end if
193    ! [2.8] Transfer pilot obs:
196    call da_read_errfac('pilot', iv % pilot_ef_u, iv % pilot_ef_v, d1, d2, d3)
198    if (iv%info(pilot)%nlocal > 0) then
199       do n = 1, iv%info(pilot)%nlocal
200          do k = 1, iv%info(pilot)%levels(n)
201             iv % pilot(n) % u(k) % error = iv % pilot(n) % u(k) % error * &
202                                            iv % pilot_ef_u
203             iv % pilot(n) % v(k) % error = iv % pilot(n) % v(k) % error * &
204                                            iv % pilot_ef_v
206          end do
207       end do
208    end if
210    ! [2.9] Transfer SSM/I obs:SSMI:
213    call da_read_errfac('ssmir', iv % ssmir_ef_speed, iv % ssmir_ef_tpw, d1, d2, d3)
215    if (iv%info(ssmi_rv)%nlocal > 0) then
216       do n = 1, iv%info(ssmi_rv)%nlocal
217          iv%ssmi_rv(n)%tpw%error = iv%ssmi_rv(n)%tpw%error * &
218                                           iv % ssmir_ef_tpw
219          iv%ssmi_rv(n)%speed%error = iv%ssmi_rv(n)%speed%error * &
220                                             iv % ssmir_ef_speed
221       end do
222    end if
225    ! iv % ssmit_ef_tb19h = 1.0 ! Tuning not yet coded.
226    ! iv % ssmit_ef_tb19v = 1.0 ! Tuning not yet coded.
227    ! iv % ssmit_ef_tb22v = 1.0 ! Tuning not yet coded.
228    ! iv % ssmit_ef_tb37h = 1.0 ! Tuning not yet coded.
229    ! iv % ssmit_ef_tb37v = 1.0 ! Tuning not yet coded.
230    ! iv % ssmit_ef_tb85h = 1.0 ! Tuning not yet coded.
231    ! iv % ssmit_ef_tb85v = 1.0 ! Tuning not yet coded.
233    if (iv%info(ssmi_tb)%nlocal > 0) then
234       ! do n = 1, iv%info(ssmi_tb)%nlocal
235       !    iv%ssmi_tb(n)%tb19h%error = iv%ssmi_tb(n)%tb19h%error
236       !    iv%ssmi_tb(n)%tb19v%error = iv%ssmi_tb(n)%tb19v%error
237       !    iv%ssmi_tb(n)%tb22v%error = iv%ssmi_tb(n)%tb22v%error
238       !    iv%ssmi_tb(n)%tb37h%error = iv%ssmi_tb(n)%tb37h%error * &
239       !                                fac_ssmit_tb37h
240       !    iv%ssmi_tb(n)%tb37v%error = iv%ssmi_tb(n)%tb37v%error * &
241       !                                fac_ssmit_tb37v
242       !    iv%ssmi_tb(n)%tb85h%error = iv%ssmi_tb(n)%tb85h%error * &
243       !                                fac_ssmit_tb85h
244       !    iv%ssmi_tb(n)%tb85v%error = iv%ssmi_tb(n)%tb85v%error * &
245       !                                fac_ssmit_tb85v
246       ! end do
247    end if
249    ! [2.10] Transfer satem obs:
251    call da_read_errfac('satem', iv % satem_ef_thickness, d1, d2, d3, d4)
253    if (iv%info(satem)%nlocal > 0) then
254       do n = 1, iv%info(satem)%nlocal
255          do k = 1, iv%info(satem)%levels(n)
256             iv % satem(n) % thickness(k) % error = iv % satem(n) % thickness(k) % error*&
257                                                    iv % satem_ef_thickness
258          end do
259       end do
260    end if
261    
262    ! [2.11] Transfer ssmt1 obs:
264    call da_read_errfac('ssmt1', iv % ssmt1_ef_t, d1, d2, d3, d4)
265       
266    if (iv%info(ssmt1)%nlocal > 0) then
267       do n = 1, iv%info(ssmt1)%nlocal
268          do k = 1, iv%info(ssmt1)%levels(n)   
269             iv % ssmt1(n) % t(k) % error = iv % ssmt1(n) % t(k) % error * &
270                                         iv % ssmt1_ef_t
271          end do
272       end do
273    end if
275    ! [2.12] Transfer ssmt2 obs:
277    call da_read_errfac('ssmt2', iv % ssmt2_ef_rh, d1, d2, d3, d4)
279    if (iv%info(ssmt2)%nlocal > 0) then
280       do n = 1, iv%info(ssmt2)%nlocal
281          do k = 1, iv%info(ssmt2)%levels(n)      
282             iv % ssmt2(n) % rh(k) % error = iv % ssmt2(n) % rh(k) % error * &
283                                          iv % ssmt2_ef_rh
284          end do
285       end do
286    end if
287    
288    ! [2.13] Transfer scatterometer obs:
290    call da_read_errfac('qscat', iv % qscat_ef_u, &
291                         iv % qscat_ef_v, d1, d2, d3)
292                            
293    if (iv%info(qscat)%nlocal > 0) then
294       do n = 1, iv%info(qscat)%nlocal
295          iv % qscat(n) % u % error = iv % qscat(n) % u % error * iv % qscat_ef_u
296          iv % qscat(n) % v % error = iv % qscat(n) % v % error * iv % qscat_ef_v
297       end do
298    end if
300    ! [2.14] Transfer profiler obs:
302    call da_read_errfac('profi', iv % profiler_ef_u, iv % profiler_ef_v, d1, d2, d3)
304    if (iv%info(profiler)%nlocal > 0) then
305       do n = 1, iv%info(profiler)%nlocal
306          do k = 1, iv%info(profiler)%levels(n)
307             iv % profiler(n) % u(k) % error = iv % profiler(n) % u(k) % error * &
308                                            iv % profiler_ef_u
309             iv % profiler(n) % v(k) % error = iv % profiler(n) % v(k) % error * &
310                                            iv % profiler_ef_v
312          end do
313       end do
314    end if
316    ! [2.15] Transfer buoy obs:
318    call da_read_errfac('buoy ', iv % buoy_ef_u, &
319                         iv % buoy_ef_v, iv % buoy_ef_t, &
320                         iv % buoy_ef_p, iv % buoy_ef_q)
321                            
322    if (iv%info(buoy)%nlocal > 0) then  
323       do n = 1, iv%info(buoy)%nlocal
324          iv % buoy(n) % u % error = iv % buoy(n) % u % error * iv % buoy_ef_u
325          iv % buoy(n) % v % error = iv % buoy(n) % v % error * iv % buoy_ef_v
326          iv % buoy(n) % t % error = iv % buoy(n) % t % error * iv % buoy_ef_t
327          iv % buoy(n) % p % error = iv % buoy(n) % p % error * iv % buoy_ef_p
328          iv % buoy(n) % q % error = iv % buoy(n) % q % error * iv % buoy_ef_q
329       end do
330    end if
332    ! [2.16] Transfer TC bogus obs:
334    call da_read_errfac('bogus', iv % bogus_ef_u, iv % bogus_ef_v, &
335                         iv % bogus_ef_t, iv % bogus_ef_q, iv % bogus_ef_slp)
337    if (iv%info(bogus)%nlocal > 0) then
338       do n = 1, iv%info(bogus)%nlocal
339          do k = 1, iv%info(bogus)%levels(n)
340             iv % bogus(n) % u(k) % error = iv % bogus(n) % u(k) % error * &
341                                            iv % bogus_ef_u
342             iv % bogus(n) % v(k) % error = iv % bogus(n) % v(k) % error * &
343                                            iv % bogus_ef_v
344             iv % bogus(n) % t(k) % error = iv % bogus(n) % t(k) % error * &
345                                            iv % bogus_ef_t
346             iv % bogus(n) % q(k) % error = iv % bogus(n) % q(k) % error * &
347                                            iv % bogus_ef_q
349          end do
351          iv % bogus(n) % slp % error = iv % bogus(n) % slp % error * iv % bogus_ef_slp
352       end do
353    end if
355    ! Transfer AIRS retrievals:
357    call da_read_errfac('airsr', iv % airsr_ef_t, iv % airsr_ef_q, d1, d3, d3)
359    if (iv%info(airsr)%nlocal > 0) then
360       do n = 1, iv%info(airsr)%nlocal
361          do k = 1, iv%info(airsr)%levels(n)
362             iv % airsr(n) % t(k) % error = iv % airsr(n) % t(k) % error * &
363                                            iv % airsr_ef_t
364             iv % airsr(n) % q(k) % error = iv % airsr(n) % q(k) % error * &
365                                            iv % airsr_ef_q
366          end do
367       end do
368    end if
370    ! Transfer mtgirs obs:
372    call da_read_errfac('mtgirs', iv % mtgirs_ef_u, iv % mtgirs_ef_v, &
373                         iv % mtgirs_ef_t, iv % mtgirs_ef_q, d1)
375    if (iv%info(mtgirs)%nlocal > 0) then
376       do n = 1, iv%info(mtgirs)%nlocal
377          do k = 1, iv%info(mtgirs)%levels(n)
378             iv % mtgirs(n) % u(k) % error = iv % mtgirs(n) % u(k) % error * &
379                                            iv % mtgirs_ef_u
380             iv % mtgirs(n) % v(k) % error = iv % mtgirs(n) % v(k) % error * &
381                                            iv % mtgirs_ef_v
382             iv % mtgirs(n) % t(k) % error = iv % mtgirs(n) % t(k) % error * &
383                                            iv % mtgirs_ef_t
384             iv % mtgirs(n) % q(k) % error = iv % mtgirs(n) % q(k) % error * &
385                                            iv % mtgirs_ef_q
387          end do
389       end do
390    end if
392    ! Transfer tamdar obs:
394    if (iv%info(tamdar)%nlocal > 0) then
395       call da_read_errfac('tamdar', iv % tamdar_ef_u, iv % tamdar_ef_v, &
396                            iv % tamdar_ef_t, iv % tamdar_ef_q, d1)
398       do n = 1, iv%info(tamdar)%nlocal
399          do k = 1, iv%info(tamdar)%levels(n)
400             iv % tamdar(n) % u(k) % error = iv % tamdar(n) % u(k) % error * &
401                                            iv % tamdar_ef_u
402             iv % tamdar(n) % v(k) % error = iv % tamdar(n) % v(k) % error * &
403                                            iv % tamdar_ef_v
404             iv % tamdar(n) % t(k) % error = iv % tamdar(n) % t(k) % error * &
405                                            iv % tamdar_ef_t
406             iv % tamdar(n) % q(k) % error = iv % tamdar(n) % q(k) % error * &
407                                            iv % tamdar_ef_q
408          end do
409       end do
410     end if
412    ! Transfer tamdar_sfc obs:
414    if (iv%info(tamdar_sfc)%nlocal > 0) then
415       do n = 1, iv%info(tamdar_sfc)%nlocal
416          iv % tamdar_sfc(n) % u % error = iv % tamdar_sfc(n) % u % error * iv % tamdar_sfc_ef_u
417          iv % tamdar_sfc(n) % v % error = iv % tamdar_sfc(n) % v % error * iv % tamdar_sfc_ef_v
418          iv % tamdar_sfc(n) % t % error = iv % tamdar_sfc(n) % t % error * iv % tamdar_sfc_ef_t
419          iv % tamdar_sfc(n) % p % error = iv % tamdar_sfc(n) % p % error * iv % tamdar_sfc_ef_p
420          iv % tamdar_sfc(n) % q % error = iv % tamdar_sfc(n) % q % error * iv % tamdar_sfc_ef_q
421       end do
422    end if
425    if (trace_use) call da_trace_exit("da_use_obs_errfac")
427 end subroutine da_use_obs_errfac